42 "This matrix is not square, the matrix MUST be square!",
43 OOMPH_CURRENT_FUNCTION,
44 OOMPH_EXCEPTION_LOCATION);
47 if (rhs.size() !=
nrow())
49 std::ostringstream error_message_stream;
50 error_message_stream <<
"The rhs vector is not the right size. It is "
51 << rhs.size() <<
", it should be " <<
nrow()
55 OOMPH_CURRENT_FUNCTION,
56 OOMPH_EXCEPTION_LOCATION);
73 Vector<std::complex<double>>& soln)
129 "This matrix is not square, the matrix MUST be square!",
130 OOMPH_CURRENT_FUNCTION,
131 OOMPH_EXCEPTION_LOCATION);
136 const double small_number = 1.0e-20;
156 for (
unsigned long i = 0;
i <
N;
i++)
159 for (
unsigned long j = 0; j <
M; j++)
161 double tmp = std::abs((*
this)(
i, j));
162 if (tmp > big) big = tmp;
167 "Singular Matrix", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
170 scaling[
i] = 1.0 / big;
187 for (
unsigned long i = 0;
i < (
N *
M);
i++)
200 for (
unsigned long j = 0; j <
M; j++)
203 unsigned long imax = 0;
205 for (
unsigned long i = 0;
i < j;
i++)
208 for (
unsigned long k = 0; k <
i; k++)
216 double largest_entry = 0.0;
217 for (
unsigned long i = j;
i <
N;
i++)
220 for (
unsigned long k = 0; k < j; k++)
226 double tmp = scaling[
i] * std::abs(sum);
227 if (tmp >= largest_entry)
237 for (
unsigned long k = 0; k <
N; k++)
239 std::complex<double> tmp =
LU_factors[
M * imax + k];
244 signature = -signature;
246 scaling[imax] = scaling[j];
258 std::complex<double> tmp = 1.0 /
LU_factors[
M * j + j];
259 for (
unsigned long i = j + 1;
i <
N;
i++)
271 std::complex<double> determinant_mantissa(1.0, 0.0);
272 std::complex<int> determinant_exponent(0, 0);
273 for (
unsigned i = 0;
i <
N;
i++)
276 std::complex<double>
entry;
277 int iexp_real = 0, iexp_imag = 0;
279 std::complex<double>(frexp(
LU_factors[
M *
i +
i].real(), &iexp_real),
285 double temp_mantissa[4];
289 temp_mantissa[0] = determinant_mantissa.real() *
entry.real();
290 temp_exp[0] = determinant_exponent.real() + iexp_real;
292 temp_mantissa[1] = determinant_mantissa.imag() *
entry.imag();
293 temp_exp[1] = determinant_exponent.imag() + iexp_imag;
296 temp_mantissa[2] = determinant_mantissa.real() *
entry.imag();
297 temp_exp[2] = determinant_exponent.real() + iexp_imag;
299 temp_mantissa[3] = determinant_mantissa.imag() *
entry.real();
300 temp_exp[3] = determinant_exponent.imag() + iexp_real;
304 for (
unsigned i = 0;
i < 3;
i += 2)
307 if (temp_exp[
i] < temp_exp[
i + 1])
310 int lower = temp_exp[
i];
313 int upper = temp_exp[
i + 1];
314 for (
int index = lower; index < upper; ++index)
316 temp_mantissa[
i] /= 2.0;
324 int lower = temp_exp[
i + 1];
327 int upper = temp_exp[
i];
328 for (
int index = lower; index < upper; ++index)
330 temp_mantissa[
i + 1] /= 2.0;
338 determinant_mantissa =
339 std::complex<double>(temp_mantissa[0] - temp_mantissa[1],
340 temp_mantissa[2] + temp_mantissa[3]);
342 determinant_exponent = std::complex<int>(temp_exp[0], temp_exp[2]);
345 determinant_mantissa =
346 std::complex<double>(frexp(determinant_mantissa.real(), &iexp_real),
347 frexp(determinant_mantissa.imag(), &iexp_imag));
349 int temp_real = determinant_exponent.real() + iexp_real;
350 int temp_imag = determinant_exponent.imag() + iexp_imag;
352 determinant_exponent = std::complex<int>(temp_real, temp_imag);
358 if (std::abs(determinant_mantissa) > 0.0)
362 if (std::abs(determinant_mantissa) < 0.0)
384 "This matrix is not square, the matrix MUST be square!",
385 OOMPH_CURRENT_FUNCTION,
386 OOMPH_EXCEPTION_LOCATION);
391 std::ostringstream error_message_stream;
392 error_message_stream <<
"The rhs vector is not the right size. It is "
393 << rhs.size() <<
", it should be " <<
N << std::endl;
396 OOMPH_CURRENT_FUNCTION,
397 OOMPH_EXCEPTION_LOCATION);
402 "Index vector has not been allocated. Have you called ludecompse()\n",
403 OOMPH_CURRENT_FUNCTION,
404 OOMPH_EXCEPTION_LOCATION);
408 std::ostringstream error_message_stream;
409 error_message_stream <<
"The Index vector is not the right size. It is "
410 <<
Index->size() <<
", it should be " <<
N
414 OOMPH_CURRENT_FUNCTION,
415 OOMPH_EXCEPTION_LOCATION);
419 unsigned long ii = 0;
420 for (
unsigned i = 0;
i <
N;
i++)
422 unsigned long ip = (*Index)[
i];
423 std::complex<double> sum = rhs[ip];
428 for (
unsigned long j = ii - 1; j <
i; j++)
441 for (
long i =
long(
N) - 1;
i >= 0;
i--)
443 std::complex<double> sum = rhs[
i];
444 for (
long j =
i + 1; j < long(
N); j++)
454 const Vector<std::complex<double>>& rhs,
455 Vector<std::complex<double>>& residual)
462 "This matrix is not square, the matrix MUST be square!",
463 OOMPH_CURRENT_FUNCTION,
464 OOMPH_EXCEPTION_LOCATION);
469 std::ostringstream error_message_stream;
470 error_message_stream <<
"The rhs vector is not the right size. It is "
471 << rhs.size() <<
", it should be " <<
N << std::endl;
474 OOMPH_CURRENT_FUNCTION,
475 OOMPH_EXCEPTION_LOCATION);
480 std::ostringstream error_message_stream;
481 error_message_stream <<
"The x vector is not the right size. It is "
482 << x.size() <<
", it should be " <<
N << std::endl;
485 OOMPH_CURRENT_FUNCTION,
486 OOMPH_EXCEPTION_LOCATION);
496 for (
unsigned long i = 0;
i <
N;
i++)
499 for (
unsigned long j = 0; j <
M; j++)
511 Vector<std::complex<double>>& soln)
517 std::ostringstream error_message_stream;
518 error_message_stream <<
"The x vector is not the right size. It is "
519 << x.size() <<
", it should be " <<
M << std::endl;
522 OOMPH_CURRENT_FUNCTION,
523 OOMPH_EXCEPTION_LOCATION);
527 if (soln.size() !=
N)
534 for (
unsigned long i = 0;
i <
N;
i++)
537 for (
unsigned long j = 0; j <
M; j++)
563 const unsigned long n = this->
nrow();
564 const unsigned long m = matrix_in.
ncol();
565 unsigned long local_nnz = 0;
569 std::complex<double>*
value = 0;
573 const int* matrix_in_row_start = matrix_in.
row_start();
574 const int* matrix_in_column_index = matrix_in.
column_index();
575 const std::complex<double>* matrix_in_value = matrix_in.
value();
578 const std::complex<double>* i_value = this->
value();
579 const int* i_row_start = this->
row_start();
590 row_start =
new int[n + 1];
594 std::set<unsigned> columns;
598 for (
unsigned long i_row = 0; i_row < n; i_row++)
601 int i_row_end = i_row_start[i_row + 1];
602 for (
int i_non_zero = i_row_start[i_row]; i_non_zero < i_row_end;
606 int matrix_in_row = i_column_index[i_non_zero];
609 int matrix_in_row_end = matrix_in_row_start[matrix_in_row + 1];
610 for (
int matrix_in_index = matrix_in_row_start[matrix_in_row];
611 matrix_in_index < matrix_in_row_end;
616 columns.insert(matrix_in_column_index[matrix_in_index]);
630 value =
new std::complex<double>[local_nnz];
634 for (
unsigned long i = 0;
i < local_nnz;
i++)
640 for (
unsigned long i_row = 0; i_row < n; i_row++)
643 int i_row_end = i_row_start[i_row + 1];
644 for (
int i_non_zero = i_row_start[i_row]; i_non_zero < i_row_end;
648 std::complex<double> i_val = i_value[i_non_zero];
651 int matrix_in_row = i_column_index[i_non_zero];
654 int matrix_in_row_end = matrix_in_row_start[matrix_in_row + 1];
655 for (
int matrix_in_index = matrix_in_row_start[matrix_in_row];
656 matrix_in_index < matrix_in_row_end;
660 int col = matrix_in_column_index[matrix_in_index];
664 for (
int insert_position =
row_start[i_row];
665 insert_position <= row_end;
668 if (insert_position == row_end)
672 std::ostringstream error_message;
673 error_message <<
"Error inserting value in result";
676 OOMPH_CURRENT_FUNCTION,
677 OOMPH_EXCEPTION_LOCATION);
683 value[insert_position] =
684 i_val * matrix_in_value[matrix_in_index];
690 value[insert_position] +=
691 i_val * matrix_in_value[matrix_in_index];
705 std::map<int, std::complex<double>>* result_maps =
706 new std::map<int, std::complex<double>>[n];
709 for (
unsigned long i_row = 0; i_row < n; i_row++)
712 int i_row_end = i_row_start[i_row + 1];
713 for (
int i_non_zero = i_row_start[i_row]; i_non_zero < i_row_end;
717 std::complex<double> i_val = i_value[i_non_zero];
720 int matrix_in_row = i_column_index[i_non_zero];
723 int matrix_in_row_end = matrix_in_row_start[matrix_in_row + 1];
724 for (
int matrix_in_index = matrix_in_row_start[matrix_in_row];
725 matrix_in_index < matrix_in_row_end;
729 int col = matrix_in_column_index[matrix_in_index];
732 result_maps[i_row][col] += i_val * matrix_in_value[matrix_in_index];
742 for (
unsigned long row = 0; row < n; row++)
744 int size = result_maps[row].size();
752 value =
new std::complex<double>[local_nnz];
756 for (
unsigned long row = 0; row < n; row++)
758 unsigned insert_position =
row_start[row];
759 for (std::map<
int, std::complex<double>>::iterator
i =
760 result_maps[row].begin();
761 i != result_maps[row].end();
765 value[insert_position] =
i->second;
770 delete[] result_maps;
778 std::vector<std::vector<int>> result_cols(n);
779 std::vector<std::vector<std::complex<double>>> result_vals(n);
782 for (
unsigned long i_row = 0; i_row < n; i_row++)
785 int i_row_end = i_row_start[i_row + 1];
786 for (
int i_non_zero = i_row_start[i_row]; i_non_zero < i_row_end;
790 std::complex<double> i_val = i_value[i_non_zero];
793 int matrix_in_row = i_column_index[i_non_zero];
796 int matrix_in_row_end = matrix_in_row_start[matrix_in_row + 1];
797 for (
int matrix_in_index = matrix_in_row_start[matrix_in_row];
798 matrix_in_index < matrix_in_row_end;
802 int col = matrix_in_column_index[matrix_in_index];
805 int size = result_cols[i_row].size();
806 for (
int i = 0;
i <= size;
i++)
811 result_cols[i_row].push_back(col);
812 result_vals[i_row].push_back(i_val *
813 matrix_in_value[matrix_in_index]);
815 else if (col == result_cols[i_row][
i])
818 result_vals[i_row][
i] +=
819 i_val * matrix_in_value[matrix_in_index];
832 for (
unsigned long row = 0; row < n; row++)
834 int size = result_cols[row].size();
842 value =
new std::complex<double>[local_nnz];
846 for (
unsigned long row = 0; row < n; row++)
848 unsigned insert_position =
row_start[row];
849 unsigned nnn = result_cols[row].size();
850 for (
unsigned i = 0;
i < nnn;
i++)
853 value[insert_position] = result_vals[row][
i];
860 const unsigned long n_nz = this->
nnz();
872 const unsigned long i_nrow = this->
nrow();
873 const unsigned long matrix_in_nrow = matrix_in.
nrow();
874 if (i_nrow != matrix_in_nrow)
876 std::ostringstream error_message;
877 error_message <<
"matrix_in has a different number of rows than\n"
880 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
883 const unsigned long i_ncol = this->
ncol();
884 const unsigned long matrix_in_ncol = matrix_in.
ncol();
885 if (i_ncol != matrix_in_ncol)
887 std::ostringstream error_message;
888 error_message <<
"matrix_in has a different number of columns than\n"
891 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
901 unsigned nrow_local = this->
nrow();
905 res_row_start.reserve(nrow_local + 1);
909 const int* i_row_start = this->
row_start();
910 const int* in_column_indices = matrix_in.
column_index();
911 const int* in_row_start = matrix_in.
row_start();
914 const std::complex<double>* i_values = this->
value();
915 const std::complex<double>* in_values = matrix_in.
value();
918 res_row_start.push_back(0);
922 for (
unsigned row_i = 0; row_i < nrow_local; row_i++)
925 std::map<int, std::complex<double>> res_row_map;
928 int i_row_end = i_row_start[row_i + 1];
929 for (
int i = i_row_start[row_i];
i < i_row_end;
i++)
931 res_row_map[i_column_indices[
i]] = i_values[
i];
935 int in_row_end = in_row_start[row_i + 1];
936 for (
int i = in_row_start[row_i];
i < in_row_end;
i++)
938 res_row_map[in_column_indices[
i]] += in_values[
i];
942 res_row_start.push_back(res_row_start.back() + res_row_map.size());
945 for (std::map<
int, std::complex<double>>::iterator it =
947 it != res_row_map.end();
950 res_column_indices.push_back(it->first);
951 res_values.push_back(it->second);
957 result_matrix.
build(res_values,
972 if (!matrix_in.
built())
974 std::ostringstream error_message;
975 error_message <<
"The matrix matrix_in is not built.\n"
976 <<
"Please build the matrix!\n";
978 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
982 const unsigned long n_row = this->
nrow();
983 const unsigned long matrix_in_nrow = matrix_in.
nrow();
984 if (n_row != matrix_in_nrow)
986 std::ostringstream error_message;
987 error_message <<
"matrix_in has a different number of rows than\n"
990 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
993 const unsigned long i_ncol = this->
ncol();
994 const unsigned long matrix_in_ncol = matrix_in.
ncol();
995 if (i_ncol != matrix_in_ncol)
997 std::ostringstream error_message;
998 error_message <<
"matrix_in has a different number of columns than\n"
1001 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1011 unsigned nrow_local = this->
nrow();
1015 res_row_start.reserve(nrow_local + 1);
1019 const int* i_row_start = this->
row_start();
1020 const int* in_column_indices = matrix_in.
column_index();
1021 const int* in_row_start = matrix_in.
row_start();
1024 const std::complex<double>* i_values = this->
value();
1025 const double* in_values = matrix_in.
value();
1028 res_row_start.push_back(0);
1032 for (
unsigned row_i = 0; row_i < nrow_local; row_i++)
1035 std::map<int, std::complex<double>> res_row_map;
1038 int i_row_end = i_row_start[row_i + 1];
1039 for (
int i = i_row_start[row_i];
i < i_row_end;
i++)
1041 res_row_map[i_column_indices[
i]] = i_values[
i];
1045 int in_row_end = in_row_start[row_i + 1];
1046 for (
int i = in_row_start[row_i];
i < in_row_end;
i++)
1048 res_row_map[in_column_indices[
i]] += in_values[
i];
1052 res_row_start.push_back(res_row_start.back() + res_row_map.size());
1055 for (std::map<
int, std::complex<double>>::iterator it =
1056 res_row_map.begin();
1057 it != res_row_map.end();
1060 res_column_indices.push_back(it->first);
1061 res_values.push_back(it->second);
1067 result_matrix.
build(res_values,
1079 const Vector<std::complex<double>>& x,
Vector<std::complex<double>>& soln)
1085 std::ostringstream error_message_stream;
1086 error_message_stream <<
"The x vector is not the right size. It is "
1087 << x.size() <<
", it should be " <<
N << std::endl;
1090 OOMPH_CURRENT_FUNCTION,
1091 OOMPH_EXCEPTION_LOCATION);
1095 if (soln.size() !=
M)
1102 for (
unsigned long i = 0;
i <
M;
i++)
1109 for (
unsigned long i = 0;
i <
N;
i++)
1111 for (
unsigned long j = 0; j <
M; j++)
1132 std::complex<double>*,
1135 std::complex<double>*,
1152 std::ostringstream error_message_stream;
1153 error_message_stream <<
"Can only solve for quadratic matrices\n"
1154 <<
"N, M " <<
N <<
" " <<
M << std::endl;
1157 OOMPH_CURRENT_FUNCTION,
1158 OOMPH_EXCEPTION_LOCATION);
1172 int nnz_aux = (int)
Nnz;
1204 if (
N != rhs.size())
1206 std::ostringstream error_message_stream;
1207 error_message_stream <<
"Size of RHS doesn't match matrix size\n"
1208 <<
"N, rhs.size() " <<
N <<
" " << rhs.size()
1212 OOMPH_CURRENT_FUNCTION,
1213 OOMPH_EXCEPTION_LOCATION);
1217 std::ostringstream error_message_stream;
1218 error_message_stream <<
"Can only solve for quadratic matrices\n"
1219 <<
"N, M " <<
N <<
" " <<
M << std::endl;
1222 OOMPH_CURRENT_FUNCTION,
1223 OOMPH_EXCEPTION_LOCATION);
1228 std::complex<double>* b =
new std::complex<double>[
N];
1231 for (
unsigned long i = 0;
i <
N;
i++)
1250 int nnz_aux = (int)
Nnz;
1270 for (
unsigned long i = 0;
i <
N;
i++)
1299 int nnz_aux = (int)
Nnz;
1325 const Vector<std::complex<double>>& rhs,
1326 Vector<std::complex<double>>& residual)
1333 "This matrix is not square, the matrix MUST be square!",
1334 OOMPH_CURRENT_FUNCTION,
1335 OOMPH_EXCEPTION_LOCATION);
1338 if (rhs.size() !=
N)
1340 std::ostringstream error_message_stream;
1341 error_message_stream <<
"The rhs vector is not the right size. It is "
1342 << rhs.size() <<
", it should be " <<
N << std::endl;
1345 OOMPH_CURRENT_FUNCTION,
1346 OOMPH_EXCEPTION_LOCATION);
1351 std::ostringstream error_message_stream;
1352 error_message_stream <<
"The x vector is not the right size. It is "
1353 << x.size() <<
", it should be " <<
N << std::endl;
1356 OOMPH_CURRENT_FUNCTION,
1357 OOMPH_EXCEPTION_LOCATION);
1361 const unsigned long r_n =
residual.size();
1368 for (
unsigned i = 0;
i <
N;
i++)
1373 for (
unsigned long j = 0; j <
N; j++)
1379 std::complex<double> a_ij =
Value[k];
1389 Vector<std::complex<double>>& soln)
1395 std::ostringstream error_message_stream;
1396 error_message_stream <<
"The x vector is not the right size. It is "
1397 << x.size() <<
", it should be " <<
M << std::endl;
1400 OOMPH_CURRENT_FUNCTION,
1401 OOMPH_EXCEPTION_LOCATION);
1405 if (soln.size() !=
N)
1410 for (
unsigned i = 0;
i <
N;
i++)
1415 for (
unsigned long j = 0; j <
N; j++)
1421 std::complex<double> a_ij =
Value[k];
1422 soln[
i] += a_ij * x[j];
1432 const Vector<std::complex<double>>& x,
Vector<std::complex<double>>& soln)
1438 std::ostringstream error_message_stream;
1439 error_message_stream <<
"The x vector is not the right size. It is "
1440 << x.size() <<
", it should be " <<
N << std::endl;
1443 OOMPH_CURRENT_FUNCTION,
1444 OOMPH_EXCEPTION_LOCATION);
1448 if (soln.size() !=
M)
1455 for (
unsigned long i = 0;
i <
M;
i++)
1461 for (
unsigned long i = 0;
i <
N;
i++)
1467 std::complex<double> a_ij =
Value[k];
1468 soln[j] += a_ij * x[
i];
1487 std::ostringstream error_message_stream;
1488 error_message_stream <<
"Can only solve for quadratic matrices\n"
1489 <<
"N, M " <<
N <<
" " <<
M << std::endl;
1492 OOMPH_CURRENT_FUNCTION,
1493 OOMPH_EXCEPTION_LOCATION);
1507 int nnz_aux = int(
Nnz);
1538 if (
N != rhs.size())
1540 std::ostringstream error_message_stream;
1541 error_message_stream <<
"Size of RHS doesn't match matrix size\n"
1542 <<
"N, rhs.size() " <<
N <<
" " << rhs.size()
1546 OOMPH_CURRENT_FUNCTION,
1547 OOMPH_EXCEPTION_LOCATION);
1551 std::ostringstream error_message_stream;
1552 error_message_stream <<
"Can only solve for quadratic matrices\n"
1553 <<
"N, M " <<
N <<
" " <<
M << std::endl;
1556 OOMPH_CURRENT_FUNCTION,
1557 OOMPH_EXCEPTION_LOCATION);
1562 std::complex<double>* b =
new std::complex<double>[
N];
1565 for (
unsigned long i = 0;
i <
N;
i++)
1583 int nnz_aux = int(
Nnz);
1603 for (
unsigned long i = 0;
i <
N;
i++)
1631 int nnz_aux = int(
Nnz);
1657 const Vector<std::complex<double>>& rhs,
1658 Vector<std::complex<double>>& residual)
1662 if (rhs.size() !=
N)
1664 std::ostringstream error_message_stream;
1665 error_message_stream <<
"The rhs vector is not the right size. It is "
1666 << rhs.size() <<
", it should be " <<
N << std::endl;
1669 OOMPH_CURRENT_FUNCTION,
1670 OOMPH_EXCEPTION_LOCATION);
1675 std::ostringstream error_message_stream;
1676 error_message_stream <<
"The x vector is not the right size. It is "
1677 << x.size() <<
", it should be " <<
M << std::endl;
1680 OOMPH_CURRENT_FUNCTION,
1681 OOMPH_EXCEPTION_LOCATION);
1690 for (
unsigned long i = 0;
i <
N;
i++)
1694 for (
long k =
Row_start[
i]; k < row_end; k++)
1697 std::complex<double> a_ij =
Value[k];
1707 Vector<std::complex<double>>& soln)
1713 std::ostringstream error_message_stream;
1714 error_message_stream <<
"The x vector is not the right size. It is "
1715 << x.size() <<
", it should be " <<
M << std::endl;
1718 OOMPH_CURRENT_FUNCTION,
1719 OOMPH_EXCEPTION_LOCATION);
1723 if (soln.size() !=
N)
1728 for (
unsigned long i = 0;
i <
N;
i++)
1732 for (
long k =
Row_start[
i]; k < row_end; k++)
1735 std::complex<double> a_ij =
Value[k];
1736 soln[
i] += a_ij * x[j];
1746 const Vector<std::complex<double>>& x,
Vector<std::complex<double>>& soln)
1752 std::ostringstream error_message_stream;
1753 error_message_stream <<
"The x vector is not the right size. It is "
1754 << x.size() <<
", it should be " <<
N << std::endl;
1757 OOMPH_CURRENT_FUNCTION,
1758 OOMPH_EXCEPTION_LOCATION);
1762 if (soln.size() !=
M)
1769 for (
unsigned long i = 0;
i <
M;
i++)
1775 for (
unsigned long i = 0;
i <
N;
i++)
1778 for (
long k =
Row_start[
i]; k < row_end; k++)
1781 std::complex<double> a_ij =
Value[k];
1782 soln[j] += a_ij * x[
i];
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
int Info
Info flag for the SuperLU solver.
void lubksub(Vector< std::complex< double >> &rhs)
LU back solve for given RHS.
int ludecompose()
LU decomposition using SuperLU.
void clean_up_memory()
LU clean up (perhaps this should happen in the destructor)
bool Doc_stats_during_solve
Flag to indicate if stats are to be displayed during solution of linear system with SuperLU.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)
Find the residulal to x of Ax=b, ie r=b-Ax.
void * F_factors
Storage for the LU factors as required by SuperLU.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
int * Column_start
Start index for column.
int * Row_index
Row index.
A class for compressed row matrices.
unsigned long nrow() const
Return the number of rows of the matrix.
void * F_factors
Storage for the LU factors as required by SuperLU.
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
unsigned long ncol() const
Return the number of columns of the matrix.
bool Doc_stats_during_solve
Flag to indicate if stats are to be displayed during solution of linear system with SuperLU.
SerialMatrixMultiplyMethod Serial_matrix_matrix_multiply_method
A switch variable for selecting the matrix multiply method for serial (non-parallel) runs....
void add(const CRDoubleMatrix &matrix_in, CRComplexMatrix &result_matrix) const
Element-wise addition of this matrix with matrix_in.
int ludecompose()
LU decomposition using SuperLU.
void lubksub(Vector< std::complex< double >> &rhs)
LU back solve for given RHS.
int Info
Info flag for the SuperLU solver.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &b, Vector< std::complex< double >> &residual)
Find the residual to x of Ax=b, i.e. r=b-Ax.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
void clean_up_memory()
LU clean up (perhaps this should happen in the destructor)
A class for compressed row matrices. This is a distributable object.
int * column_index()
Access to C-style column index array.
int * row_start()
Access to C-style row_start array.
unsigned long ncol() const
Return the number of columns of the matrix.
bool built() const
access function to the Built flag - indicates whether the matrix has been build - i....
double * value()
Access to C-style value array.
unsigned long nrow() const
Return the number of rows of the matrix.
int * row_start()
Access to C-style row_start array.
int * Row_start
Start index for row.
int * Column_index
Column index.
void build(const Vector< T > &value, const Vector< int > &column_index, const Vector< int > &row_start, const unsigned long &n, const unsigned long &m)
Build matrix from compressed representation. Number of nonzero entries is read off from value,...
int * column_index()
Access to C-style column index array.
void build_without_copy(T *value, int *column_index, int *row_start, const unsigned long &nnz, const unsigned long &n, const unsigned long &m)
Function to build matrix from pointers to arrays which hold the row starts, column indices and non-ze...
virtual int ludecompose()
LU decomposition of the matrix using the appropriate linear solver. Return the sign of the determinan...
virtual unsigned long ncol() const =0
Return the number of columns of the matrix.
virtual void lubksub(Vector< std::complex< double >> &rhs)
LU backsubstitue a previously LU-decomposed matrix; The passed rhs will be over-written with the solu...
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.
virtual void solve(Vector< std::complex< double >> &rhs)
Complete LU solve (replaces matrix by its LU decomposition and overwrites RHS with solution)....
bool Overwrite_matrix_storage
Boolean flag used to decided whether the LU decomposition is stored separately, or not.
void lubksub(Vector< std::complex< double >> &rhs)
Overload the LU back substitution. Note that the rhs will be overwritten with the solution vector.
virtual ~DenseComplexMatrix()
Destructor, delete the storage for Index vector and LU storage (if any)
Vector< long > * Index
Pointer to storage for the index of permutations in the LU solve.
void residual(const Vector< std::complex< double >> &x, const Vector< std::complex< double >> &rhs, Vector< std::complex< double >> &residual)
Find the residual of Ax=b, ie r=b-Ax for the "solution" x.
void multiply(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the matrix by the vector x: soln=Ax.
int ludecompose()
Overload the LU decomposition. For this storage scheme the matrix storage will be over-writeen by its...
void delete_lu_factors()
Function that deletes the storage for the LU_factors, if it is not the same as the matrix storage.
std::complex< double > * LU_factors
Pointer to storage for the LU decomposition.
void multiply_transpose(const Vector< std::complex< double >> &x, Vector< std::complex< double >> &soln)
Multiply the transposed matrix by the vector x: soln=A^T x.
std::complex< double > & entry(const unsigned long &i, const unsigned long &j)
The access function that will be called by the read-write round-bracket operator.
unsigned long N
Number of rows.
std::complex< double > * Matrixdata
Internal representation of matrix as a pointer to data.
unsigned long M
Number of columns.
An OomphLibError object which should be thrown when an run-time error is encountered....
unsigned long Nnz
Number of non-zero values (i.e. size of Value array)
T * value()
Access to C-style value array.
T * Value
Internal representation of the matrix values, a pointer.
unsigned long nnz() const
Return the number of nonzero entries.
unsigned long N
Number of rows.
unsigned long M
Number of columns.
A slight extension to the standard template vector class so that we can include "graceful" array rang...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
int superlu_complex(int *, int *, int *, int *, std::complex< double > *, int *, int *, std::complex< double > *, int *, int *, int *, void *, int *)
///////////////////////////////////////////////////////////////////// ///////////////////////////////...