21 #include "libmesh/sparse_matrix.h" 24 #include "libmesh/dof_map.h" 25 #include "libmesh/dense_matrix.h" 26 #include "libmesh/diagonal_matrix.h" 27 #include "libmesh/laspack_matrix.h" 28 #include "libmesh/libmesh_logging.h" 29 #include "libmesh/eigen_sparse_matrix.h" 30 #include "libmesh/parallel.h" 31 #include "libmesh/petsc_matrix.h" 32 #include "libmesh/trilinos_epetra_matrix.h" 33 #include "libmesh/numeric_vector.h" 34 #include "libmesh/enum_solver_package.h" 38 #ifdef LIBMESH_HAVE_GZSTREAM 39 # include "gzstream.h" 46 #ifdef LIBMESH_HAVE_CXX11_REGEX 66 _use_hash_table(false)
92 const std::vector<numeric_index_type> & brows,
93 const std::vector<numeric_index_type> & bcols)
95 libmesh_assert_equal_to (dm.
m() / brows.size(), dm.
n() / bcols.size());
98 (dm.
m() / brows.size());
100 libmesh_assert_equal_to (dm.
m()%blocksize, 0);
101 libmesh_assert_equal_to (dm.
n()%blocksize, 0);
103 std::vector<numeric_index_type> rows, cols;
105 rows.reserve(blocksize*brows.size());
106 cols.reserve(blocksize*bcols.size());
108 for (
auto & row : brows)
112 for (
unsigned int v=0; v<blocksize; v++)
116 for (
auto & col : bcols)
120 for (
unsigned int v=0; v<blocksize; v++)
124 this->add_matrix (dm, rows, cols);
137 libmesh_not_implemented();
140 os <<
"Real part:" << std::endl;
144 os << std::setw(8) << (*this)(i,j).
real() <<
" ";
148 os << std::endl <<
"Imaginary part:" << std::endl;
152 os << std::setw(8) << (*this)(i,j).
imag() <<
" ";
163 template <
typename T>
164 std::unique_ptr<SparseMatrix<T>>
173 return std::make_unique<DiagonalMatrix<T>>(comm);
176 switch (solver_package)
179 #ifdef LIBMESH_HAVE_LASPACK 181 return std::make_unique<LaspackMatrix<T>>(comm);
185 #ifdef LIBMESH_HAVE_PETSC 187 return std::make_unique<PetscMatrix<T>>(comm);
191 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA 193 return std::make_unique<EpetraMatrix<T>>(comm);
197 #ifdef LIBMESH_HAVE_EIGEN 199 return std::make_unique<EigenSparseMatrix<T>>(comm);
203 libmesh_error_msg(
"ERROR: Unrecognized solver package: " << solver_package);
208 template <
typename T>
213 this->vector_mult_add(dest,arg);
218 template <
typename T>
229 template <
typename T>
233 libmesh_not_implemented();
237 template <
typename T>
246 template <
typename T>
249 parallel_object_only();
254 end_dof = this->row_stop();
258 if (this->processor_id() == 0)
260 libmesh_assert_equal_to (first_dof, 0);
268 if (c != static_cast<T>(0.0))
270 os << i <<
" " << j <<
" " << c << std::endl;
277 os << (*
this)(i,j) <<
" ";
282 std::vector<numeric_index_type> ibuf, jbuf;
287 this->comm().receive(p, ibuf);
288 this->comm().receive(p, jbuf);
289 this->comm().receive(p, cbuf);
290 libmesh_assert_equal_to (ibuf.size(), jbuf.size());
291 libmesh_assert_equal_to (ibuf.size(), cbuf.size());
295 libmesh_assert_greater_equal (ibuf.front(), currenti);
296 libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
298 std::size_t currentb = 0;
299 for (;currenti <= ibuf.back(); ++currenti)
305 if (currentb < ibuf.size() &&
306 ibuf[currentb] == currenti &&
309 os << currenti <<
" " << j <<
" " << cbuf[currentb] << std::endl;
318 if (currentb < ibuf.size() &&
319 ibuf[currentb] == currenti &&
322 os << cbuf[currentb] <<
" ";
326 os << static_cast<T>(0.0) <<
" ";
334 for (; currenti != this->m(); ++currenti)
337 os << static_cast<T>(0.0) <<
" ";
344 std::vector<numeric_index_type> ibuf, jbuf;
354 if (c != static_cast<T>(0.0))
362 this->comm().send(0,ibuf);
363 this->comm().send(0,jbuf);
364 this->comm().send(0,cbuf);
369 template <
typename T>
372 parallel_object_only();
377 end_dof = this->row_stop();
381 if (this->processor_id() == 0)
383 std::unique_ptr<std::ofstream> file;
386 file = std::make_unique<std::ofstream>(
name.c_str());
390 std::size_t sparsity_nonzeros = this->n_nonzeros();
392 std::size_t real_nonzeros = 0;
394 libmesh_assert_equal_to(first_dof, 0);
400 if (c != static_cast<T>(0.0))
408 std::size_t nonzeros_on_p = 0;
409 this->comm().receive(p, nonzeros_on_p);
410 real_nonzeros += nonzeros_on_p;
413 if (sparsity_nonzeros &&
414 sparsity_nonzeros != real_nonzeros)
415 libmesh_warning(sparsity_nonzeros <<
416 " nonzeros allocated, but " <<
417 real_nonzeros <<
" used.");
423 os <<
"%Mat Object: () " << this->n_processors() <<
" MPI processes\n" 424 <<
"% type: " << (this->n_processors() > 1 ?
"mpi" :
"seq") <<
"aij\n" 425 <<
"% Size = " << this->m() <<
' ' << this->n() <<
'\n' 426 <<
"% Nonzeros = " << real_nonzeros <<
'\n' 427 <<
"zzz = zeros(" << real_nonzeros <<
",3);\n" 437 if (c != static_cast<T>(0.0))
440 os << (i+1) <<
' ' << (j+1) <<
" " << c <<
'\n';
445 std::vector<numeric_index_type> ibuf, jbuf;
449 this->comm().receive(p, ibuf);
450 this->comm().receive(p, jbuf);
451 this->comm().receive(p, cbuf);
452 libmesh_assert_equal_to (ibuf.size(), jbuf.size());
453 libmesh_assert_equal_to (ibuf.size(), cbuf.size());
456 os << ibuf[n] <<
' ' << jbuf[n] <<
" " << cbuf[n] <<
'\n';
459 os <<
"];\n" <<
"Mat_sparse = spconvert(zzz);" << std::endl;
463 std::vector<numeric_index_type> ibuf, jbuf;
465 std::size_t my_nonzeros = 0;
474 if (c != static_cast<T>(0.0))
483 this->comm().send(0,my_nonzeros);
484 this->comm().send(0,ibuf);
485 this->comm().send(0,jbuf);
486 this->comm().send(0,cbuf);
492 template <
typename T>
495 libmesh_not_implemented_msg
496 (
"libMesh cannot write PETSc binary-format files from non-PETSc matrices");
501 template <
typename T>
504 libmesh_not_implemented_msg
505 (
"libMesh cannot write PETSc HDF5-format files from non-PETSc matrices");
510 template <
typename T>
514 std::ifstream in (filename.c_str());
516 (!in.good(),
"ERROR: cannot read file:\n\t" <<
522 const bool gzipped_file = (basename.rfind(
".gz") == basename.size() - 3);
525 basename.remove_suffix(3);
527 if (basename.rfind(
".matlab") == basename.size() - 7 ||
528 basename.rfind(
".m") == basename.size() - 2)
529 this->read_matlab(filename);
530 else if (basename.rfind(
".petsc64") == basename.size() - 8)
532 #ifndef LIBMESH_HAVE_PETSC 533 libmesh_error_msg(
"Cannot load PETSc matrix file " <<
534 filename <<
" without PETSc-enabled libMesh.");
536 #if LIBMESH_DOF_ID_BYTES != 8 537 libmesh_error_msg(
"Cannot load 64-bit PETSc matrix file " <<
538 filename <<
" with non-64-bit libMesh.");
540 this->read_petsc_binary(filename);
542 else if (basename.rfind(
".petsc32") == basename.size() - 8)
544 #ifndef LIBMESH_HAVE_PETSC 545 libmesh_error_msg(
"Cannot load PETSc matrix file " <<
546 filename <<
" without PETSc-enabled libMesh.");
548 #if LIBMESH_DOF_ID_BYTES != 4 549 libmesh_error_msg(
"Cannot load 32-bit PETSc matrix file " <<
550 filename <<
" with non-32-bit libMesh.");
552 this->read_petsc_binary(filename);
555 libmesh_error_msg(
" ERROR: Unrecognized matrix file extension on: " 557 <<
"\n I understand the following:\n\n" 558 <<
" *.matlab -- Matlab sparse matrix format\n" 559 <<
" *.m -- Matlab sparse matrix format\n" 560 <<
" *.petsc32 -- PETSc binary format, 32-bit\n" 561 <<
" *.petsc64 -- PETSc binary format, 64-bit\n" 566 template <
typename T>
569 LOG_SCOPE(
"read_matlab()",
"SparseMatrix");
571 #ifndef LIBMESH_HAVE_CXX11_REGEX 572 libmesh_not_implemented();
575 parallel_object_only();
577 const bool gzipped_file = (filename.rfind(
".gz") == filename.size() - 3);
585 std::vector<numeric_index_type> new_row_starts, new_row_stops,
586 new_col_starts, new_col_stops;
589 new_col_start, new_col_stop;
600 std::unique_ptr<std::istream> file;
609 std::vector<std::tuple<numeric_index_type, numeric_index_type, T>> entries;
615 if (this->processor_id() == 0)
619 const std::regex start_regex
620 (
"\\s*\\w+\\s*=\\s*\\[");
621 const std::regex end_regex
626 #ifdef LIBMESH_HAVE_GZSTREAM 627 auto inf = std::make_unique<igzstream>();
629 inf->open(filename.c_str(), std::ios::in);
630 file = std::move(inf);
632 libmesh_error_msg(
"ERROR: need gzstream to handle .gz files!!!");
637 auto inf = std::make_unique<std::ifstream>();
642 inf->open(new_name.c_str(), std::ios::in);
643 file = std::move(inf);
648 const std::regex size_regex
649 (
"%\\s*[Ss][Ii][Zz][Ee]\\s*=\\s*(\\d+)\\s+(\\d+)");
650 const std::string whitespace =
" \t";
652 bool have_started =
false;
653 bool have_ended =
false;
654 std::size_t largest_i_seen = 0, largest_j_seen = 0;
658 std::size_t current_row = 1;
660 for (std::string line; std::getline(*file, line);)
670 std::istringstream l(line);
675 l >> i >> j >>
value;
680 (!have_started,
"Confused by premature entries in matrix file " << filename);
682 entries.emplace_back(cast_int<numeric_index_type>(i),
683 cast_int<numeric_index_type>(j),
687 (!i || !j,
"Expected 1-based indexing in matrix file " 690 current_row = std::max(current_row, i);
694 "Can't handle out-of-order entries in matrix file " 697 largest_i_seen = std::max(i, largest_i_seen);
698 largest_j_seen = std::max(j, largest_j_seen);
701 else if (std::regex_search(line, sm, size_regex))
703 const std::string msize = sm[1];
704 const std::string nsize = sm[2];
705 m = std::stoull(msize);
706 n = std::stoull(nsize);
709 else if (std::regex_search(line, start_regex))
712 else if (std::regex_search(line, end_regex))
720 (!have_started,
"Confused by missing assignment beginning in matrix file " << filename);
723 (!have_ended,
"Confused by missing assignment ending in matrix file " << filename);
726 (m > largest_i_seen,
"Confused by missing final row(s) in matrix file " << filename);
729 (m > 0 && m < largest_i_seen,
"Confused by extra final row(s) in matrix file " << filename);
735 (n > largest_j_seen,
"Confused by missing final column(s) in matrix file " << filename);
738 (n > 0 && n < largest_j_seen,
"Confused by extra final column(s) in matrix file " << filename);
743 this->comm().broadcast(m);
744 this->comm().broadcast(n);
748 this->comm().broadcast(m);
749 this->comm().broadcast(n);
756 new_row_start = this->row_start(),
757 new_row_stop = this->row_stop();
759 new_col_start = this->col_start(),
760 new_col_stop = this->col_stop();
765 new_row_start = this->processor_id() * m / this->n_processors(),
766 new_row_stop = (this->processor_id()+1) * m / this->n_processors();
768 new_col_start = this->processor_id() * n / this->n_processors(),
769 new_col_stop = (this->processor_id()+1) * n / this->n_processors();
772 this->comm().gather(0, new_row_start, new_row_starts);
773 this->comm().gather(0, new_row_stop, new_row_stops);
774 this->comm().gather(0, new_col_start, new_col_starts);
775 this->comm().gather(0, new_col_stop, new_col_stops);
784 off_diagonal_nonzeros =0;
786 if (this->processor_id() == 0)
795 for (
auto [i, j,
value] : entries)
801 while (current_row >= (new_row_stops[current_proc]+1))
803 current_on_diagonal_nonzeros = 0;
804 current_off_diagonal_nonzeros = 0;
808 if (j >= (new_col_starts[current_proc]+1) &&
809 j < (new_col_stops[current_proc]+1))
811 ++current_on_diagonal_nonzeros;
812 on_diagonal_nonzeros =
813 std::max(on_diagonal_nonzeros,
814 current_on_diagonal_nonzeros);
818 ++current_off_diagonal_nonzeros;
819 off_diagonal_nonzeros =
820 std::max(off_diagonal_nonzeros,
821 current_off_diagonal_nonzeros);
826 this->comm().broadcast(on_diagonal_nonzeros);
827 this->comm().broadcast(off_diagonal_nonzeros);
830 new_row_stop-new_row_start,
831 new_col_stop-new_col_start,
832 on_diagonal_nonzeros,
833 off_diagonal_nonzeros);
838 if (this->processor_id() == 0)
839 for (
auto [i, j,
value] : entries)
840 this->
set(i-1, j-1,
value);
848 template <
typename T>
851 libmesh_not_implemented_msg
852 (
"libMesh cannot read PETSc binary-format files into non-PETSc matrices");
857 template <
typename T>
860 libmesh_not_implemented_msg
861 (
"libMesh cannot read PETSc HDF5-format files into non-PETSc matrices");
866 template <
typename T>
871 for (
const auto i :
make_range(this->row_start(), this->row_stop()))
872 for (
const auto j :
make_range(this->col_start(), this->col_stop()))
873 this->
set(i, j, (*this)(i, j) *
scale);
878 template <
typename T>
881 auto diff_mat = this->clone();
882 diff_mat->add(-1.0, other_mat);
883 return diff_mat->l1_norm();
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
bool closed()
Checks that the library has been closed.
virtual std::size_t n_nonzeros() const
SparseMatrix(const Parallel::Communicator &comm)
Constructor; initializes the matrix to be empty, without any structure, i.e.
This helper class can be called on multiple threads to compute the sparsity pattern (or graph) of the...
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
The libMesh namespace provides an interface to certain functionality in the library.
const SparsityPattern::Build * get_sparsity_pattern() const
uint8_t processor_id_type
virtual void zero()=0
Set all entries to zero.
This class handles the numbering of degrees of freedom on a mesh.
void libmesh_ignore(const Args &...)
dof_id_type numeric_index_type
bool _is_initialized
Flag that tells if init() has been called.
auto l1_norm_diff(const NumericVector< T > &vec1, const NumericVector< T > &vec2)
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
std::string unzip_file(std::string_view name)
Create an unzipped copy of a bz2 or xz file, returning the name of the now-unzipped file that can be ...
An object whose state is distributed along a set of processors.
std::string_view basename_of(const std::string &fullname)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
bool initialized()
Checks that library initialization has been done.
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
SolverPackage
Defines an enum for various linear solver packages.
MatrixBuildType
Defines an enum for matrix build types.
Defines a dense matrix for use in Finite Element-type computations.
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...