62#pragma GCC diagnostic push
63#pragma GCC diagnostic ignored "-Wold-style-cast"
64#pragma GCC diagnostic ignored "-Wignored-attributes"
68#include <unsupported/Eigen/CXX11/Tensor>
69#pragma GCC diagnostic pop
71#define MAXBUFSIZE (static_cast<int> (1e6))
72#define PBSTR "||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||"
82#pragma GCC diagnostic push
83#pragma GCC diagnostic ignored "-Wold-style-cast"
85#pragma GCC diagnostic pop
100template<
typename Type>
115template<
typename T,
int dim>
116void exportMatrix(Eigen::Matrix < T, -1, dim > & matrix, word Name,
117 word type =
"python", word folder =
"./Matrices");
130void exportMatrix(List <Eigen::MatrixXd>& matrix, word Name,
131 word type =
"python", word folder =
"./Matrices");
145void exportTensor(Eigen::Tensor<T, 3 > tensor, word Name,
146 word type =
"python", word folder =
"./Matrices");
165template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
166void read_fields(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
168 fileName casename,
int first_snap = 0,
int n_snap = 0);
188template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
189void read_fields(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
191 GeometricField<Type, PatchField, GeoMesh>& field,
192 fileName casename,
int first_snap = 0,
int n_snap = 0);
207template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
210 GeometricField<Type, PatchField, GeoMesh>& field,
226template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
228 PtrList<GeometricField<Type, PatchField, GeoMesh>>& Lfield,
229 GeometricField<Type, PatchField, GeoMesh>& field,
245template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
248 const GeometricField<Type, PatchField, GeoMesh>& field,
249 const fileName casename);
264template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
265void readLastFields(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
267 const GeometricField<Type, PatchField, GeoMesh>& field,
268 const fileName casename);
281template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
282void exportFields(PtrList<GeometricField<Type, PatchField, GeoMesh>>&
284 word folder, word fieldname);
308List <Eigen::MatrixXd>
readMatrix(word folder, word mat_name);
322template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
324 fileName subfolder, fileName folder,
338template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
340 fileName subfolder, fileName folder);
352void exportList(
T& list, word folder, word filename);
366template <
typename T,
int Nrows,
typename IND>
368 word folder, word MatrixName);
381template <
typename T,
int Nrows,
typename IND>
383 word folder, word MatrixName);
394template <
typename MatrixType>
406template <
typename MatrixType>
419template <
typename TensorType>
432template <
typename TensorType>
444template <
typename MatrixType>
459template <
typename T,
int Nrows,
typename IND>
461 word folder, word MatrixName);
472template <
typename MatrixType>
486template <
typename MatrixType>
498int numberOfFiles(word folder, word MatrixName, word ext =
"");
507void writePoints(pointField points, fileName folder, fileName subfolder);
517void save(
const List<Eigen::SparseMatrix<T>>& MatrixList, word folder,
521void load(List<Eigen::SparseMatrix<T>>& MatrixList, word folder,
529Ostream&
operator<< (Ostream& os,
const Eigen::Matrix < T, -1, -1 > & mat)
531 os << mat.rows() << mat.cols() << UList<T>(
const_cast<Eigen::MatrixXd&
>
532 (mat).data(), mat.size());
537Istream&
operator>> (Istream& is, Eigen::Matrix < T, -1, -1 > & mat)
541 mat.resize(nrow, ncol);
542 UList<T> list(mat.data(), nrow * ncol);
548Ostream&
operator<< (Ostream& os,
const Eigen::Matrix < T, -1, 1 > & mat)
550 os << mat.rows() << mat.cols() << UList<T>(
const_cast<Eigen::VectorXd&
>
551 (mat).data(), mat.size());
556Istream&
operator>> (Istream& is, Eigen::Matrix < T, -1, 1 > & mat)
560 mat.resize(nrow, ncol);
561 UList<T> list(mat.data(), nrow * ncol);
567Ostream&
operator<< (Ostream& os,
const Eigen::Tensor<T, 3>& tens)
569 os << tens.dimension(0) << tens.dimension(1) << tens.dimension(
570 2) << UList<T>(
const_cast<Eigen::Tensor<T, 3>&
>(tens).data(),
579 is >> d1 >> d2 >> d3;
580 tens.resize(d1, d2, d3);
581 UList<T> list(tens.data(), d1 * d2 * d3);
588template <
typename T,
int whatever,
typename IND>
590 word folder, word MatrixName)
592 typedef Eigen::Triplet<int> Trip;
593 std::vector<Trip> res;
596 std::fstream writeFile;
597 writeFile.open(folder + MatrixName, std::ios::binary | std::ios::out);
599 if (writeFile.is_open())
601 IND rows, cols, nnzs, outS, innS;
604 nnzs = m.nonZeros() ;
605 outS = m.outerSize();
606 innS = m.innerSize();
607 writeFile.write(
reinterpret_cast<const char*
> (&rows),
sizeof(IND));
608 writeFile.write(
reinterpret_cast<const char*
> (&cols),
sizeof(IND));
609 writeFile.write(
reinterpret_cast<const char*
> (&nnzs),
sizeof(IND));
610 writeFile.write(
reinterpret_cast<const char*
> (&outS),
sizeof(IND));
611 writeFile.write(
reinterpret_cast<const char*
> (&innS),
sizeof(IND));
612 writeFile.write(
reinterpret_cast<const char*
>(m.valuePtr()),
613 sizeof(
T ) * m.nonZeros());
614 writeFile.write(
reinterpret_cast<const char*
>(m.outerIndexPtr()),
615 sizeof(IND) * m.outerSize());
616 writeFile.write(
reinterpret_cast<const char*
>(m.innerIndexPtr()),
617 sizeof(IND) * m.nonZeros());
622template <
typename T,
int whatever,
typename IND>
624 word folder, word MatrixName)
626 std::fstream readFile;
627 readFile.open(folder + MatrixName, std::ios::binary | std::ios::in);
629 if (!readFile.good())
631 std::cout << folder + MatrixName <<
632 " file does not exist, try to rerun the Offline Stage!" << std::endl;
635 else if (readFile.is_open())
637 IND rows, cols, nnz, inSz, outSz;
638 readFile.read(
reinterpret_cast<char*
>(&rows ),
sizeof(IND));
639 readFile.read(
reinterpret_cast<char*
>(&cols ),
sizeof(IND));
640 readFile.read(
reinterpret_cast<char*
>(&nnz ),
sizeof(IND));
641 readFile.read(
reinterpret_cast<char*
>(&inSz ),
sizeof(IND));
642 readFile.read(
reinterpret_cast<char*
>(&outSz),
sizeof(IND));
643 m.resize(rows, cols);
645 m.resizeNonZeros(nnz);
646 readFile.read(
reinterpret_cast<char*
>(m.valuePtr()),
sizeof(
T ) * nnz );
647 readFile.read(
reinterpret_cast<char*
>(m.outerIndexPtr()),
sizeof(IND) * outSz);
648 readFile.read(
reinterpret_cast<char*
>(m.innerIndexPtr()),
sizeof(IND) * nnz );
654template <
typename MatrixType>
659 std::ofstream out(folder + MatrixName,
660 std::ios::out | std::ios::binary | std::ios::trunc);
661 typename MatrixType::Index rows = Matrix.rows(), cols = Matrix.cols();
662 out.write(
reinterpret_cast<char*
> (&rows),
sizeof(
typename MatrixType::Index));
663 out.write(
reinterpret_cast<char*
> (&cols),
sizeof(
typename MatrixType::Index));
664 out.write(
reinterpret_cast<char*
> (Matrix.data()),
665 rows * cols *
sizeof(
typename MatrixType::Scalar) );
669template <
typename MatrixType>
674 in.open(folder + MatrixName, std::ios::in | std::ios::binary);
678 std::cout << folder + MatrixName <<
679 " file does not exist, try to rerun the Offline Stage!" << std::endl;
682 else if (in.is_open())
684 typename MatrixType::Index rows = 0, cols = 0;
685 in.read(
reinterpret_cast<char*
> (&rows),
sizeof(
typename MatrixType::Index));
686 in.read(
reinterpret_cast<char*
> (&cols),
sizeof(
typename MatrixType::Index));
687 Matrix.resize(rows, cols);
688 in.read(
reinterpret_cast<char*
>(Matrix.data()),
689 rows * cols *
sizeof(
typename MatrixType::Scalar) );
695template <
typename TensorType>
699 std::ofstream out(folder + MatrixName,
700 std::ios::out | std::ios::binary | std::ios::trunc);
701 typename TensorType::Dimensions dim = Tensor.dimensions();
704 for (
int k = 0; k < dim.size(); k++)
709 out.write(
reinterpret_cast<char*
> (&dim),
710 sizeof(
typename TensorType::Dimensions));
711 out.write(
reinterpret_cast<char*
> (Tensor.data()),
712 tot *
sizeof(
typename TensorType::Scalar) );
716template <
typename TensorType>
721 in.open(folder + MatrixName, std::ios::in | std::ios::binary);
722 typename TensorType::Dimensions dim;
723 in.read(
reinterpret_cast<char*
> (&dim),
724 sizeof(
typename TensorType::Dimensions));
725 auto dims = Tensor.dimensions();
727 "The rank of the tensor you want to fill does not coincide with the rank of the tensor you are reading");
730 for (
int k = 0; k < dim.size(); k++)
736 in.read(
reinterpret_cast<char*
>(Tensor.data()),
737 tot *
sizeof(
typename TensorType::Scalar) );
741template <
typename MatrixType>
743 word folder, word MatrixName)
745 for (
int i = 0;
i < MatrixList.size();
i++)
751template <
typename MatrixType>
753 word folder, word MatrixName)
755 for (
int i = 0;
i < MatrixList.size();
i++)
761template <
typename T,
int whatever,
typename IND>
763 List<Eigen::SparseMatrix<T, whatever, IND>>& m, word folder, word MatrixName)
765 int number_of_files = numberOfFiles(folder, MatrixName);
766 std::cout <<
"Reading the Matrix " + folder + MatrixName << std::endl;
768 "Check if the file you are trying to read exists" );
769 Eigen::SparseMatrix<T, whatever, IND>
A;
771 for (
int i = 0;
i < number_of_files;
i++)
773 ReadSparseMatrix(
A, folder, MatrixName + name(
i));
779template <
typename MatrixType>
784 std::cout <<
"Reading the Matrix " + folder + MatrixName << std::endl;
786 "Check if the file you are trying to read exists" );
789 for (
int i = 0;
i < number_of_files;
i++)
Header file of the EigenFunctions class.
Implementation of the assert function for ITHACA-FV.
#define M_Assert(Expr, Msg)
Header file of the ITHACAutilities namespace.
Ostream & operator<<(Ostream &os, const Eigen::Matrix< T, -1, -1 > &mat)
Istream & operator>>(Istream &is, Eigen::Matrix< T, -1, -1 > &mat)
Namespace for input-output manipulation.
void writePoints(pointField points, fileName folder, fileName subfolder)
Write points of a mesh to a file.
void ReadDenseMatrixList(List< MatrixType > &m, word folder, word MatrixName)
Read a sparse matrix list to bynary files.
void SaveSparseMatrix(Eigen::SparseMatrix< T, Nrows, IND > &m, word folder, word MatrixName)
Export an Eigen sparse matrix into bynary format file.
void ReadDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Read a dense matrix from a binary format file.
int numberOfFiles(word folder, word MatrixName, word ext)
Count the number of files with a certain name prefix, used by reading functions of list of matrices.
void readConvergedFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, GeometricField< Type, PatchField, GeoMesh > &field, fileName casename)
Function to read a list of volVectorField from name of the field including only converged snapshots.
void save(const List< Eigen::SparseMatrix< T > > &MatrixList, word folder, word MatrixName)
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
void load(List< Eigen::SparseMatrix< T > > &MatrixList, word folder, word MatrixName)
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void read_last_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, const GeometricField< Type, PatchField, GeoMesh > &field, const fileName casename)
Function to read a list of fields from the name of the field including only the last snapshot.
void readMiddleFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, GeometricField< Type, PatchField, GeoMesh > &field, fileName casename)
Funtion to read a list of volVectorField from name of the field including all the intermediate snapsh...
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
void exportList(T &list, word folder, word filename)
Export a list to file.
void ReadSparseMatrixList(List< Eigen::SparseMatrix< T, Nrows, IND > > &m, word folder, word MatrixName)
Read a sparse matrix list to bynary files.
void SaveDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Save a dense matrix to a binary format file.
void ReadDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Read a dense tensor from file.
void exportFvMatrix(fvMatrix< Type > &Matrix, word folder, word MatrixName)
Export an Fv Matrix to folder together with its source term.
void SaveDenseMatrixList(List< MatrixType > &MatrixList, word folder, word MatrixName)
Save a dense matrix list to bynary files.
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
void SaveSparseMatrixList(List< MatrixType > &MatrixList, word folder, word MatrixName)
Save a sparse matrix list to bynary files.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
void readLastFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, const GeometricField< Type, PatchField, GeoMesh > &field, const fileName casename)
Funtion to read a list of volVectorField from name of the field including only the last snapshots.
void exportTensor(Eigen::Tensor< T, 3 > tensor, word Name, word type, word folder)
Export the reduced tensor in numpy (tipo=python), matlab (tipo=matlab) and txt (tipo=eigen) format.
void ReadSparseMatrix(Eigen::SparseMatrix< T, Nrows, IND > &m, word folder, word MatrixName)
Read an Eigen sparse matrix from a bynary format file.
void SaveDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Save a dense tensor to file.
void printProgress(double percentage)
Print progress bar given the percentage.