Loading...
Searching...
No Matches
ITHACAstream.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 ITHACAstream
26Description
27 Methods for input and outputs from and to files
28SourceFiles
29 ITHACAstream.C
30 ITHACAstreamTemplates.C
31\*---------------------------------------------------------------------------*/
32
38
39#ifndef ITHACAstream_H
40#define ITHACAstream_H
41
42#include "fvCFD.H"
43#include "IOmanip.H"
44#include <stdio.h>
45#include <sys/types.h>
46#include <dirent.h>
47#include <algorithm>
48#include <fstream>
49#include <string>
50#include <stdexcept>
51#include <sstream>
52#include <vector>
53#include <cstdio>
54#include <typeinfo>
55#include <iostream>
56#include <cassert>
57#include <zlib.h>
58#include <map>
59#include "ITHACAassert.H"
60#include "ITHACAparameters.H"
61#include "ITHACAutilities.H"
62#pragma GCC diagnostic push
63#pragma GCC diagnostic ignored "-Wold-style-cast"
64#pragma GCC diagnostic ignored "-Wignored-attributes"
65#include <Eigen/Eigen>
66#include "EigenFunctions.H"
67#include "colormod.H"
68
69#include <unsupported/Eigen/CXX11/Tensor>
70#pragma GCC diagnostic pop
71
72#define MAXBUFSIZE (static_cast<int> (1e6))
73#define PBSTR "||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||"
74#define PBWIDTH 60
75
76// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77
78
79/*---------------------------------------------------------------------------*\
80 Namespace ITHACAstream Declaration
81\*---------------------------------------------------------------------------*/
82
83#pragma GCC diagnostic push
84#pragma GCC diagnostic ignored "-Wold-style-cast"
85#include "cnpy.H"
86#pragma GCC diagnostic pop
87
89namespace ITHACAstream
90{
91
92//--------------------------------------------------------------------------
101template<typename Type>
102void exportFvMatrix(fvMatrix<Type>& Matrix, word folder,
103 word MatrixName);
104
105//--------------------------------------------------------------------------
116template<typename T, int dim>
117void exportMatrix(Eigen::Matrix < T, -1, dim > & matrix, word Name,
118 word type = "python", word folder = "./Matrices");
119
120//--------------------------------------------------------------------------
131void exportMatrix(List <Eigen::MatrixXd>& matrix, word Name,
132 word type = "python", word folder = "./Matrices");
133
134//--------------------------------------------------------------------------
145template<typename T>
146void exportTensor(Eigen::Tensor<T, 3 > tensor, word Name,
147 word type = "python", word folder = "./Matrices");
148
149
150//----------------------------------------------------------------------
166template<class Type, template<class> class PatchField, class GeoMesh>
167void read_fields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
168 Lfield, word Name,
169 fileName casename, int first_snap = 0, int n_snap = 0);
170
171
172
173//----------------------------------------------------------------------
189template<class Type, template<class> class PatchField, class GeoMesh >
190void read_fields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
191 Lfield,
192 GeometricField<Type, PatchField, GeoMesh>& field,
193 fileName casename, int first_snap = 0, int n_snap = 0);
194
195//----------------------------------------------------------------------
208template<class Type, template<class> class PatchField, class GeoMesh>
209void readMiddleFields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
210 Lfield,
211 GeometricField<Type, PatchField, GeoMesh>& field,
212 fileName casename);
213
214//----------------------------------------------------------------------
227template<class Type, template<class> class PatchField, class GeoMesh>
229 PtrList<GeometricField<Type, PatchField, GeoMesh >>& Lfield,
230 GeometricField<Type, PatchField, GeoMesh>& field,
231 fileName casename);
232
233//----------------------------------------------------------------------
241
246template<class Type, template<class> class PatchField, class GeoMesh>
247void read_last_fields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
248 Lfield,
249 const GeometricField<Type, PatchField, GeoMesh>& field,
250 const fileName casename);
251
252//----------------------------------------------------------------------
265template<class Type, template<class> class PatchField, class GeoMesh>
266void readLastFields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
267 Lfield,
268 const GeometricField<Type, PatchField, GeoMesh>& field,
269 const fileName casename);
270
271//--------------------------------------------------------------------------
282template<class Type, template<class> class PatchField, class GeoMesh>
283void exportFields(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
284 field,
285 word folder, word fieldname);
286
287//--------------------------------------------------------------------------
289/* One has to provide the complete filename with the absolute or relative path */
295Eigen::MatrixXd readMatrix(word filename);
296
297//--------------------------------------------------------------------------
299/* One has to provide the folder containing the matrix files and the filename of the
300 the matrix. Since it is stored as a List of matrices each matrix must be stored in a
301 different file with the following format:
302 matFileName(i)_mat.txt */
309List <Eigen::MatrixXd> readMatrix(word folder, word mat_name);
310
311//--------------------------------------------------------------------------
323template<class Type, template<class> class PatchField, class GeoMesh>
324void exportSolution(GeometricField<Type, PatchField, GeoMesh>& s,
325 fileName subfolder, fileName folder,
326 word fieldName);
327
328//--------------------------------------------------------------------------
339template<class Type, template<class> class PatchField, class GeoMesh>
340void exportSolution(GeometricField<Type, PatchField, GeoMesh>& s,
341 fileName subfolder, fileName folder);
342
343//--------------------------------------------------------------------------
352template<typename T>
353void exportList(T& list, word folder, word filename);
354
355
356//--------------------------------------------------------------------------
367template <typename T, int Nrows, typename IND>
368void SaveSparseMatrix(Eigen::SparseMatrix<T, Nrows, IND>& m,
369 word folder, word MatrixName);
370
371//--------------------------------------------------------------------------
382template <typename T, int Nrows, typename IND>
383void ReadSparseMatrix(Eigen::SparseMatrix<T, Nrows, IND>& m,
384 word folder, word MatrixName);
385
386//--------------------------------------------------------------------------
395template <typename MatrixType>
396void SaveDenseMatrix(MatrixType& Matrix, word folder, word MatrixName);
397
398//--------------------------------------------------------------------------
407template <typename MatrixType>
408void ReadDenseMatrix(MatrixType& Matrix, word folder, word MatrixName);
409
410
411//----------------------------------------------------------------------
420template <typename TensorType>
421void SaveDenseTensor(TensorType& Tensor, word folder, word MatrixName);
422
423
424//----------------------------------------------------------------------
433template <typename TensorType>
434void ReadDenseTensor(TensorType& Tensor, word folder, word MatrixName);
435
436//--------------------------------------------------------------------------
445template <typename MatrixType>
446void SaveSparseMatrixList(List<MatrixType>& MatrixList, word folder,
447 word MatrixName);
448
449//--------------------------------------------------------------------------
460template <typename T, int Nrows, typename IND>
461void ReadSparseMatrixList(List<Eigen::SparseMatrix<T, Nrows, IND >>& m,
462 word folder, word MatrixName);
463
464//--------------------------------------------------------------------------
473template <typename MatrixType>
474void SaveDenseMatrixList(List<MatrixType>& MatrixList, word folder,
475 word MatrixName);
476
477
478//--------------------------------------------------------------------------
487template <typename MatrixType>
488void ReadDenseMatrixList(List<MatrixType>& m, word folder,
489 word MatrixName);
490
491//--------------------------------------------------------------------------
499int numberOfFiles(word folder, word MatrixName, word ext = "");
500
501//--------------------------------------------------------------------------
508void writePoints(pointField points, fileName folder, fileName subfolder);
509
510//--------------------------------------------------------------------------
515void printProgress(double percentage);
516
517template<typename T>
518void save(const List<Eigen::SparseMatrix<T >>& MatrixList, word folder,
519 word MatrixName);
520
521template<typename T>
522void load(List<Eigen::SparseMatrix<T >>& MatrixList, word folder,
523 word MatrixName);
524
525//--------------------------------------------------------------------------
538template<class Type, template<class> class PatchField, class GeoMesh >
539GeometricField<Type, PatchField, GeoMesh> readFieldByIndex(
540 const GeometricField<Type, PatchField, GeoMesh>& field,
541 fileName casename,
542 label index
543);
544
545template <typename T>
546void read_snapshot(T& snapshot, const Foam::word snap_time,
547 Foam::word path, Foam::word name = "default_name");
548
549
550};
551
552namespace Foam
553{
554template<typename T>
555Ostream& operator<< (Ostream& os, const Eigen::Matrix < T, -1, -1 > & mat)
556{
557 os << mat.rows() << mat.cols() << UList<T>(const_cast<Eigen::Matrix<T, -1, -1>&>
558 (mat).data(), mat.size());
559 return os;
560}
561
562template<typename T>
563Istream& operator>> (Istream& is, Eigen::Matrix < T, -1, -1 > & mat)
564{
565 int nrow, ncol;
566 is >> nrow >> ncol;
567 mat.resize(nrow, ncol);
568 UList<T> list(mat.data(), nrow * ncol);
569 is >> list;
570 return is;
571}
572
573template<typename T>
574Ostream& operator<< (Ostream& os, const Eigen::Matrix < T, -1, 1 > & mat)
575{
576 os << mat.rows() << mat.cols() << UList<T>(const_cast<Eigen::Matrix<T, -1, 1>&>
577 (mat).data(), mat.size());
578 return os;
579}
580
581template<typename T>
582Istream& operator>> (Istream& is, Eigen::Matrix < T, -1, 1 > & mat)
583{
584 int nrow, ncol;
585 is >> nrow >> ncol;
586 mat.resize(nrow, ncol);
587 UList<T> list(mat.data(), nrow * ncol);
588 is >> list;
589 return is;
590}
591
592// Support for Eigen matrix transpose expression
593template<typename T>
594Ostream& operator<< (Ostream& os, const Eigen::Transpose<Eigen::Matrix<T, -1, -1>>& mat)
595{
596 Eigen::Matrix<T, -1, -1> evaluated = mat;
597 os << evaluated.rows() << evaluated.cols() << UList<T>(evaluated.data(), evaluated.size());
598 return os;
599}
600
601// Support for Eigen vector transpose expression
602template<typename T>
603Ostream& operator<< (Ostream& os, const Eigen::Transpose<Eigen::Matrix<T, -1, 1>>& mat)
604{
605 Eigen::Matrix<T, 1, -1> evaluated = mat;
606 os << evaluated.rows() << evaluated.cols() << UList<T>(evaluated.data(), evaluated.size());
607 return os;
608}
609
610// Support for Eigen::Block (column/row views from matrices)
611template<typename T>
612Ostream& operator<< (Ostream& os, const Eigen::Block<Eigen::Matrix<T, -1, -1>, -1, 1, true>& mat)
613{
614 Eigen::Matrix<T, -1, 1> evaluated = mat;
615 os << evaluated;
616 return os;
617}
618
619template<typename T>
620Ostream& operator<< (Ostream& os, const Eigen::Block<Eigen::Matrix<T, -1, -1>, 1, -1, false>& mat)
621{
622 Eigen::Matrix<T, 1, -1> evaluated = mat;
623 os << evaluated;
624 return os;
625}
626
627template<typename T>
628Ostream& operator<< (Ostream& os, const Eigen::Tensor<T, 3>& tens)
629{
630 os << tens.dimension(0) << tens.dimension(1) << tens.dimension(
631 2) << UList<T>(const_cast<Eigen::Tensor<T, 3>&>(tens).data(),
632 tens.size());
633 return os;
634}
635
636template<typename T>
637Istream& operator>> (Istream& is, Eigen::Tensor<T, 3>& tens)
638{
639 int d1, d2, d3;
640 is >> d1 >> d2 >> d3;
641 tens.resize(d1, d2, d3);
642 UList<T> list(tens.data(), d1 * d2 * d3);
643 is >> list;
644 return is;
645}
646
647inline Ostream& operator<<(Ostream& os, const Color::Modifier& mod)
648{
649 // Output ANSI code to Foam Ostream
650 os << "\033[" << int(mod.getCode()) << "m";
651 return os;
652}
653
654}
655
656template <typename T, int whatever, typename IND>
657void ITHACAstream::SaveSparseMatrix(Eigen::SparseMatrix<T, whatever, IND>& m,
658 word folder, word MatrixName)
659{
660 typedef Eigen::Triplet<int> Trip;
661 std::vector<Trip> res;
662 m.makeCompressed();
663 mkDir(folder);
664 std::fstream writeFile;
665 writeFile.open(folder + MatrixName, std::ios::binary | std::ios::out);
666
667 if (writeFile.is_open())
668 {
669 IND rows, cols, nnzs, outS, innS;
670 rows = m.rows() ;
671 cols = m.cols() ;
672 nnzs = m.nonZeros() ;
673 outS = m.outerSize();
674 innS = m.innerSize();
675 writeFile.write(reinterpret_cast<const char*> (& rows), sizeof(IND));
676 writeFile.write(reinterpret_cast<const char*> (& cols), sizeof(IND));
677 writeFile.write(reinterpret_cast<const char*> (& nnzs), sizeof(IND));
678 writeFile.write(reinterpret_cast<const char*> (& outS), sizeof(IND));
679 writeFile.write(reinterpret_cast<const char*> (& innS), sizeof(IND));
680 writeFile.write(reinterpret_cast<const char*>(m.valuePtr()),
681 sizeof(T ) * m.nonZeros());
682 writeFile.write(reinterpret_cast<const char*>(m.outerIndexPtr()),
683 sizeof(IND) * m.outerSize());
684 writeFile.write(reinterpret_cast<const char*>(m.innerIndexPtr()),
685 sizeof(IND) * m.nonZeros());
686 writeFile.close();
687 }
688}
689
690template <typename T, int whatever, typename IND>
691void ITHACAstream::ReadSparseMatrix(Eigen::SparseMatrix<T, whatever, IND>& m,
692 word folder, word MatrixName)
693{
694 std::fstream readFile;
695 readFile.open(folder + MatrixName, std::ios::binary | std::ios::in);
696
697 if (!readFile.good())
698 {
699 Foam::Info << folder + MatrixName <<
700 " file does not exist, try to rerun the Offline Stage!" << Foam::endl;
701 exit(EXIT_FAILURE);
702 }
703 else if (readFile.is_open())
704 {
705 IND rows, cols, nnz, inSz, outSz;
706 readFile.read(reinterpret_cast<char*>(& rows ), sizeof(IND));
707 readFile.read(reinterpret_cast<char*>(& cols ), sizeof(IND));
708 readFile.read(reinterpret_cast<char*>(& nnz ), sizeof(IND));
709 readFile.read(reinterpret_cast<char*>(& inSz ), sizeof(IND));
710 readFile.read(reinterpret_cast<char*>(& outSz), sizeof(IND));
711 m.resize(rows, cols);
712 m.makeCompressed();
713 m.resizeNonZeros(nnz);
714 readFile.read(reinterpret_cast<char*>(m.valuePtr()), sizeof(T ) * nnz );
715 readFile.read(reinterpret_cast<char*>(m.outerIndexPtr()), sizeof(IND) * outSz);
716 readFile.read(reinterpret_cast<char*>(m.innerIndexPtr()), sizeof(IND) * nnz );
717 m.finalize();
718 readFile.close();
719 }
720}
721
722template <typename MatrixType>
723void ITHACAstream::SaveDenseMatrix(MatrixType& Matrix, word folder,
724 word MatrixName)
725{
726 mkDir(folder);
727 std::ofstream out(folder + MatrixName,
728 std::ios::out | std::ios::binary | std::ios::trunc);
729 typename MatrixType::Index rows = Matrix.rows(), cols = Matrix.cols();
730 out.write(reinterpret_cast<char*> (& rows),
731 sizeof(typename MatrixType::Index));
732 out.write(reinterpret_cast<char*> (& cols),
733 sizeof(typename MatrixType::Index));
734 out.write(reinterpret_cast<char*> (Matrix.data()),
735 rows * cols * sizeof(typename MatrixType::Scalar) );
736 out.close();
737}
738
739template <typename MatrixType>
740void ITHACAstream::ReadDenseMatrix(MatrixType& Matrix, word folder,
741 word MatrixName)
742{
743 std::ifstream in;
744 in.open(folder + MatrixName, std::ios::in | std::ios::binary);
745
746 if (!in.good())
747 {
748 Foam::Info << folder + MatrixName <<
749 " file does not exist, try to rerun the Offline Stage!" << Foam::endl;
750 exit(EXIT_FAILURE);
751 }
752 else if (in.is_open())
753 {
754 typename MatrixType::Index rows = 0, cols = 0;
755 in.read(reinterpret_cast<char*> (& rows), sizeof(typename MatrixType::Index));
756 in.read(reinterpret_cast<char*> (& cols), sizeof(typename MatrixType::Index));
757 Matrix.resize(rows, cols);
758 in.read( reinterpret_cast<char*>(Matrix.data()),
759 rows * cols * sizeof(typename MatrixType::Scalar) );
760 in.close();
761 }
762}
763
764
765template <typename TensorType>
766void ITHACAstream::SaveDenseTensor(TensorType& Tensor, word folder,
767 word MatrixName)
768{
769 std::ofstream out(folder + MatrixName,
770 std::ios::out | std::ios::binary | std::ios::trunc);
771 typename TensorType::Dimensions dim = Tensor.dimensions();
772 int tot = 1;
773
774 for (int k = 0; k < dim.size(); k++)
775 {
776 tot *= dim[k];
777 }
778
779 out.write(reinterpret_cast<char*> (& dim),
780 sizeof(typename TensorType::Dimensions));
781 out.write(reinterpret_cast<char*> (Tensor.data()),
782 tot * sizeof(typename TensorType::Scalar) );
783 out.close();
784}
785
786template <typename TensorType>
787void ITHACAstream::ReadDenseTensor(TensorType& Tensor, word folder,
788 word MatrixName)
789{
790 std::ifstream in;
791 in.open(folder + MatrixName, std::ios::in | std::ios::binary);
792 typename TensorType::Dimensions dim;
793 in.read(reinterpret_cast<char*> (& dim),
794 sizeof(typename TensorType::Dimensions));
795 auto dims = Tensor.dimensions();
796 M_Assert(dims.size() == dim.size(),
797 "The rank of the tensor you want to fill does not coincide with the rank of the tensor you are reading");
798 int tot = 1;
799
800 for (int k = 0; k < dim.size(); k++)
801 {
802 tot *= dim[k];
803 }
804
805 Tensor.resize(dim);
806 in.read( reinterpret_cast<char*>(Tensor.data()),
807 tot * sizeof(typename TensorType::Scalar) );
808 in.close();
809}
810
811template <typename MatrixType>
812void ITHACAstream::SaveSparseMatrixList(List<MatrixType>& MatrixList,
813 word folder, word MatrixName)
814{
815 for (int i = 0; i < MatrixList.size(); i++)
816 {
817 SaveSparseMatrix(MatrixList[i], folder, MatrixName + name(i));
818 }
819}
820
821template <typename MatrixType>
822void ITHACAstream::SaveDenseMatrixList(List<MatrixType>& MatrixList,
823 word folder, word MatrixName)
824{
825 for (int i = 0; i < MatrixList.size(); i++)
826 {
827 SaveDenseMatrix(MatrixList[i], folder, MatrixName + name(i));
828 }
829}
830
831template <typename T, int whatever, typename IND>
833 List<Eigen::SparseMatrix<T, whatever, IND >>& m, word folder, word MatrixName)
834{
835 int number_of_files = numberOfFiles(folder, MatrixName);
836 Foam::Info << "Reading the Matrix " + folder + MatrixName << Foam::endl;
837 M_Assert(number_of_files != 0,
838 "Check if the file you are trying to read exists" );
839 Eigen::SparseMatrix<T, whatever, IND> A;
840
841 for (int i = 0; i < number_of_files; i++)
842 {
843 ReadSparseMatrix(A, folder, MatrixName + name(i));
844 m.append(A);
845 }
846}
847
848
849template <typename MatrixType>
850void ITHACAstream::ReadDenseMatrixList(List<MatrixType>& m, word folder,
851 word MatrixName)
852{
853 int number_of_files = numberOfFiles(folder, MatrixName);
854 Foam::Info << "Reading the Matrix " + folder + MatrixName << Foam::endl;
855 M_Assert(number_of_files != 0,
856 "Check if the file you are trying to read exists" );
857 MatrixType A;
858
859 for (int i = 0; i < number_of_files; i++)
860 {
861 ReadDenseMatrix(A, folder, MatrixName + name(i));
862 m.append(A);
863 }
864}
865
866
867
868
869
870// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
871
872
873
874#endif
875
876
877
878
879
880
Header file of the EigenFunctions class.
Implementation of the assert function for ITHACA-FV.
Header file of the ITHACAutilities namespace.
Simple header and source file of the Color::Modifier class to change color to the output stream.
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.
GeometricField< Type, PatchField, GeoMesh > readFieldByIndex(const GeometricField< Type, PatchField, GeoMesh > &field, fileName casename, label index)
Function to read a single field by index from a folder.
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 exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
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.