Loading...
Searching...
No Matches
Foam2Eigen.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 -------------------------------------------------------------------------------
12 License
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/>.
24 Class
25 foam2eigen
26 Description
27 Converts foam objects into Eigen and viceversa
28 SourceFiles
29 foam2eigen.C
30 \*---------------------------------------------------------------------------*/
31
36
37#ifndef Foam2Eigen_H
38#define Foam2Eigen_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "pointMesh.H"
43#include "pointPatchField.H"
44#include "ITHACAassert.H"
45#include "ITHACAutilities.H"
46#include <tuple>
47#include <sys/stat.h>
48#pragma GCC diagnostic push
49#pragma GCC diagnostic ignored "-Wold-style-cast"
50#include <Eigen/Eigen>
51#include <type_traits>
52#pragma GCC diagnostic pop
53// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54
55
56/*---------------------------------------------------------------------------*\
57 Class foam2eigen Declaration
58\*---------------------------------------------------------------------------*/
59
62{
63 private:
64
65 public:
66
67 //----------------------------------------------------------------------
89 template <class type_foam_matrix, class type_A, class type_B>
90 static void fvMatrix2Eigen(fvMatrix<type_foam_matrix> foam_matrix, type_A& A,
91 type_B& b);
92 /*
93 template <typename ScalarType, typename SparseMatType, typename VecType>
94 static void fvMat2Eigen(fvMatrix<ScalarType> foam_matrix,
95 SparseMatType& A,
96 VecType& b);*/
97 template<typename SparseMatType, typename VecType>
98 static void fvMat2Eigen(Foam::fvMatrix<scalar> foam_matrix,
99 SparseMatType& A,
100 VecType& b);
101
102 template<typename SparseMatType, typename VecType>
103 static void fvMat2Eigen(Foam::fvMatrix<vector> foam_matrix,
104 SparseMatType& A,
105 VecType& b);
106
107 //----------------------------------------------------------------------
123 template <class type_foam_matrix, class type_A>
124 static void fvMatrix2EigenM(fvMatrix<type_foam_matrix>& foam_matrix,
125 type_A& A);
126
127 //----------------------------------------------------------------------
141 template <class type_foam_matrix, class type_B>
142 static void fvMatrix2EigenV(fvMatrix<type_foam_matrix>& foam_matrix,
143 type_B& b);
144
145
146 //----------------------------------------------------------------------
160 template<class Type, template<class> class PatchField, class GeoMesh>
161 static Eigen::MatrixXd PtrList2Eigen(
162 PtrList<GeometricField<Type, PatchField, GeoMesh >>& fields,
163 label Nfields = -1);
164
165 //----------------------------------------------------------------------
175 template<class Type, template<class> class PatchField, class GeoMesh>
176 static Eigen::MatrixXd field2Eigen(
177 GeometricField<Type, PatchField, GeoMesh>& field);
178
179
180 //----------------------------------------------------------------------
190 template<class Type, template<class> class PatchField, class GeoMesh>
191 static Eigen::Map<Eigen::MatrixXd> field2EigenMap(
192 GeometricField<Type, PatchField, GeoMesh>& field);
193
194 //----------------------------------------------------------------------
205 template<class Type, template<class> class PatchField, class GeoMesh>
206 static Eigen::Map<Eigen::MatrixXd> field2EigenMapBC(
207 GeometricField<Type, PatchField, GeoMesh>&
208 field, int BC_index);
209
210 //----------------------------------------------------------------------
220 template<class Type, class GeoMesh>
221 static Eigen::VectorXd field2Eigen(const
222 DimensionedField<Type, GeoMesh>& field);
223
224 //----------------------------------------------------------------------
234 template<class Type>
235 static Eigen::VectorXd field2Eigen(const Field<Type>& field);
236
237 //----------------------------------------------------------------------
248 template<template<class> class PatchField, class GeoMesh>
249 static List<Eigen::VectorXd> field2EigenBC(
250 GeometricField<scalar, PatchField, GeoMesh>& field);
251
252 //----------------------------------------------------------------------
263 template<template<class> class PatchField, class GeoMesh>
264 static List<Eigen::VectorXd> field2EigenBC(
265 GeometricField<vector, PatchField, GeoMesh>& field);
266
267 //----------------------------------------------------------------------
278 template<template<class> class PatchField, class GeoMesh>
279 static List<Eigen::VectorXd> field2EigenBC(
280 GeometricField<tensor, PatchField, GeoMesh>& field);
281
282
283 //----------------------------------------------------------------------
295 template<template<class> class PatchField, class GeoMesh>
296 static List<Eigen::MatrixXd> PtrList2EigenBC(
297 PtrList<GeometricField<scalar, PatchField, GeoMesh >>& fields,
298 label Nfields = -1);
299
300 //----------------------------------------------------------------------
312 template<template<class> class PatchField, class GeoMesh >
313 static List<Eigen::MatrixXd> PtrList2EigenBC(
314 PtrList<GeometricField<vector, PatchField, GeoMesh >>& fields,
315 label Nfields = -1);
316
317 //----------------------------------------------------------------------
329 template<template<class> class PatchField, class GeoMesh >
330 static List<Eigen::MatrixXd> PtrList2EigenBC(
331 PtrList<GeometricField<tensor, PatchField, GeoMesh >>& fields,
332 label Nfields = -1);
333
334 //----------------------------------------------------------------------
346 template<template<class> class PatchField, class GeoMesh >
347 static GeometricField<scalar, PatchField, GeoMesh> Eigen2field(
348 GeometricField<scalar, PatchField, GeoMesh>& field,
349 Eigen::VectorXd& eigen_vector, bool correctBC = true);
350
351 //----------------------------------------------------------------------
363 template<template<class> class PatchField, class GeoMesh>
364 static GeometricField<vector, PatchField, GeoMesh> Eigen2field(
365 GeometricField<vector, PatchField, GeoMesh>& field,
366 Eigen::VectorXd& eigen_vector, bool correctBC = true);
367
368 //----------------------------------------------------------------------
381 template<template<class> class PatchField, class GeoMesh>
382 static GeometricField<scalar, PatchField, GeoMesh> Eigen2field(
383 GeometricField<scalar, PatchField, GeoMesh>& field,
384 Eigen::VectorXd& eigen_vector, List<Eigen::VectorXd>& eigen_vector_boundary);
385
386 //----------------------------------------------------------------------
399 template<template<class> class PatchField, class GeoMesh>
400 static GeometricField<vector, PatchField, GeoMesh> Eigen2field(
401 GeometricField<vector, PatchField, GeoMesh>& field,
402 Eigen::VectorXd& eigen_vector, List<Eigen::VectorXd>& eigen_vector_boundary);
403
404
405 //----------------------------------------------------------------------
417 template<template<class> class PatchField, class GeoMesh>
418 static GeometricField<tensor, PatchField, GeoMesh> Eigen2field(
419 GeometricField<tensor, PatchField, GeoMesh>& field,
420 Eigen::VectorXd& eigen_vector, bool correctBC = true);
421
422 //----------------------------------------------------------------------
433 template <class Type>
434 static Field<Type> Eigen2field(
435 Field<Type>& field, Eigen::MatrixXd& matrix, bool correctBC = true);
436
437 //----------------------------------------------------------------------
450 template<class Type, template<class> class PatchField, class GeoMesh>
451 static std::tuple<Eigen::MatrixXd, Eigen::VectorXd> projectFvMatrix(
452 fvMatrix<Type>& matrix,
453 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes, label Nmodes);
454
455
456 //----------------------------------------------------------------------
482 template<class Type, template<class> class PatchField, class GeoMesh >
483 static Eigen::VectorXd projectField(
484 GeometricField<Type, PatchField, GeoMesh>& field,
485 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
486 label Nmodes);
487
488 //----------------------------------------------------------------------
500 template<class Type, template<class> class PatchField, class GeoMesh >
501 static Eigen::MatrixXd MassMatrix(
502 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes, label Nmodes);
503
504 //----------------------------------------------------------------------
516 template<class Type>
517 static std::tuple<List<Eigen::SparseMatrix<double >>, List<Eigen::VectorXd >>
518 LFvMatrix2LSM(PtrList<fvMatrix<Type >>& MatrixList);
519
520 //--------------------------------------------------------------------------
529 template <class type_matrix>
530 static Eigen::Matrix<type_matrix, Eigen::Dynamic, Eigen::Dynamic>
531 List2EigenMatrix ( List<type_matrix> list );
532
533 //--------------------------------------------------------------------------
542 template <class type_matrix>
543 static List<type_matrix> EigenMatrix2List (
544 Eigen::Matrix<type_matrix, Eigen::Dynamic, Eigen::Dynamic> matrix );
545
546 //----------------------------------------------------------------------
555 template <class type_list>
556 static Eigen::MatrixXd field2Eigen(const List<type_list>& list);
557};
558
559#endif
Implementation of the assert function for ITHACA-FV.
Header file of the ITHACAutilities namespace.
Class to converts OpenFOAM objects into Eigen and viceversa.
Definition Foam2Eigen.H:62
static GeometricField< scalar, PatchField, GeoMesh > Eigen2field(GeometricField< scalar, PatchField, GeoMesh > &field, Eigen::VectorXd &eigen_vector, bool correctBC=true)
Convert a vector in Eigen format into an OpenFOAM scalar GeometricField.
Definition Foam2Eigen.C:533
static List< Eigen::VectorXd > field2EigenBC(GeometricField< scalar, PatchField, GeoMesh > &field)
Convert an OpenFOAM scalar field to a List of Eigen Vectors, one for each boundary.
Definition Foam2Eigen.C:253
static Field< Type > Eigen2field(Field< Type > &field, Eigen::MatrixXd &matrix, bool correctBC=true)
Converts a matrix in Eigen format into an OpenFOAM Field.
static std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > LFvMatrix2LSM(PtrList< fvMatrix< Type > > &MatrixList)
Convert a PtrList of OpenFOAM fvMatrix into a tuple of lists of Eigen Sparse Matrices and source vect...
static Eigen::Matrix< type_matrix, Eigen::Dynamic, Eigen::Dynamic > List2EigenMatrix(List< type_matrix > list)
Convert a Foam List into an Eigen matrix with one column.
static Eigen::VectorXd projectField(GeometricField< Type, PatchField, GeoMesh > &field, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes)
Perform the projection of an OpenFOAM field onto a set of modes using the Eigen matrix multiplication...
static Eigen::MatrixXd MassMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes)
Obtain the Mass Matrix from a list of reduced basis.
static List< type_matrix > EigenMatrix2List(Eigen::Matrix< type_matrix, Eigen::Dynamic, Eigen::Dynamic > matrix)
Convert an Eigen matrix with one column into a Foam List.
static Eigen::Map< Eigen::MatrixXd > field2EigenMapBC(GeometricField< Type, PatchField, GeoMesh > &field, int BC_index)
Map a scalar OpenFOAM field boundary into an Eigen Matrix.
static void fvMatrix2EigenV(fvMatrix< type_foam_matrix > &foam_matrix, type_B &b)
Convert a ldu OpenFOAM matrix into a source vector b.
static Eigen::MatrixXd field2Eigen(GeometricField< Type, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
static List< Eigen::MatrixXd > PtrList2EigenBC(PtrList< GeometricField< scalar, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert an OpenFOAM scalar field to a List of Eigen Vectors, one for each boundary.
Definition Foam2Eigen.C:273
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field).
Definition Foam2Eigen.C:665
static Eigen::Map< Eigen::MatrixXd > field2EigenMap(GeometricField< Type, PatchField, GeoMesh > &field)
Convert a scalar OpenFOAM field into an Eigen Vector.
static void fvMatrix2Eigen(fvMatrix< type_foam_matrix > foam_matrix, type_A &A, type_B &b)
Convert a FvMatrix OpenFOAM matrix (Linear System) into a Eigen Matrix A and a source vector b.
static std::tuple< Eigen::MatrixXd, Eigen::VectorXd > projectFvMatrix(fvMatrix< Type > &matrix, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes)
Obtain the Mass Matrix from a list of reduced basis.
static Eigen::VectorXd field2Eigen(const DimensionedField< Type, GeoMesh > &field)
Convert an OpenFOAM field into an Eigen Vector.
static Eigen::VectorXd field2Eigen(const Field< Type > &field)
Convert an OpenFOAM field into an Eigen Vector.
static Eigen::MatrixXd field2Eigen(const List< type_list > &list)
Function to convert an OpenFOAM list to an Eigen Matrix.
static void fvMatrix2EigenM(fvMatrix< type_foam_matrix > &foam_matrix, type_A &A)
Convert a ldu OpenFOAM matrix into a Eigen Matrix A.