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 "ITHACAassert.H"
43#include "ITHACAutilities.H"
44#include <tuple>
45#include <sys/stat.h>
46#pragma GCC diagnostic push
47#pragma GCC diagnostic ignored "-Wold-style-cast"
48#include <Eigen/Eigen>
49#pragma GCC diagnostic pop
50// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52
53/*---------------------------------------------------------------------------*\
54 Class foam2eigen Declaration
55\*---------------------------------------------------------------------------*/
56
59{
60 private:
61
62 public:
63
64 //----------------------------------------------------------------------
86 template <class type_foam_matrix, class type_A, class type_B>
87 static void fvMatrix2Eigen(fvMatrix<type_foam_matrix> foam_matrix, type_A& A,
88 type_B& b);
89
90 //----------------------------------------------------------------------
106 template <class type_foam_matrix, class type_A>
107 static void fvMatrix2EigenM(fvMatrix<type_foam_matrix>& foam_matrix,
108 type_A& A);
109
110 //----------------------------------------------------------------------
124 template <class type_foam_matrix, class type_B>
125 static void fvMatrix2EigenV(fvMatrix<type_foam_matrix>& foam_matrix,
126 type_B& b);
127
128
129 //----------------------------------------------------------------------
143 template<class Type, template<class> class PatchField, class GeoMesh>
144 static Eigen::MatrixXd PtrList2Eigen(
145 PtrList<GeometricField<Type, PatchField, GeoMesh >> & fields,
146 label Nfields = -1);
147
148 //----------------------------------------------------------------------
158 template<template<class> class PatchField, class GeoMesh >
159 static Eigen::VectorXd field2Eigen(GeometricField<tensor, PatchField, GeoMesh>&
160 field);
161
162 //----------------------------------------------------------------------
172 template<template<class> class PatchField, class GeoMesh>
173 static Eigen::VectorXd field2Eigen(GeometricField<vector, PatchField, GeoMesh>&
174 field);
175
176 //----------------------------------------------------------------------
186 template<template<class> class PatchField, class GeoMesh>
187 static Eigen::VectorXd field2Eigen(GeometricField<scalar, PatchField, GeoMesh>&
188 field);
189
190 //----------------------------------------------------------------------
200 template<template<class> class PatchField, class GeoMesh>
201 static Eigen::Map<Eigen::MatrixXd> field2EigenMap(
202 GeometricField<scalar, PatchField, GeoMesh>&
203 field);
204
205 //----------------------------------------------------------------------
216 template<template<class> class PatchField, class GeoMesh>
217 static Eigen::Map<Eigen::MatrixXd> field2EigenMapBC(
218 GeometricField<scalar, PatchField, GeoMesh>&
219 field, int BC_index);
220
221 //----------------------------------------------------------------------
231 template<class Type, class GeoMesh>
232 static Eigen::VectorXd field2Eigen(const
233 DimensionedField<Type, GeoMesh>& field);
234
235 //----------------------------------------------------------------------
245 template<class Type>
246 static Eigen::VectorXd field2Eigen(const Field<Type>& field);
247
248 //----------------------------------------------------------------------
259 template<template<class> class PatchField, class GeoMesh>
260 static List<Eigen::VectorXd> field2EigenBC(
261 GeometricField<scalar, PatchField, GeoMesh>& field);
262
263 //----------------------------------------------------------------------
274 template<template<class> class PatchField, class GeoMesh>
275 static List<Eigen::VectorXd> field2EigenBC(
276 GeometricField<vector, PatchField, GeoMesh>& field);
277
278 //----------------------------------------------------------------------
289 template<template<class> class PatchField, class GeoMesh>
290 static List<Eigen::VectorXd> field2EigenBC(
291 GeometricField<tensor, PatchField, GeoMesh>& field);
292
293
294 //----------------------------------------------------------------------
306 template<template<class> class PatchField, class GeoMesh>
307 static List<Eigen::MatrixXd> PtrList2EigenBC(
308 PtrList<GeometricField<scalar, PatchField, GeoMesh >> & fields,
309 label Nfields = -1);
310
311 //----------------------------------------------------------------------
323 template<template<class> class PatchField, class GeoMesh >
324 static List<Eigen::MatrixXd> PtrList2EigenBC(
325 PtrList<GeometricField<vector, PatchField, GeoMesh >> & fields,
326 label Nfields = -1);
327
328 //----------------------------------------------------------------------
340 template<template<class> class PatchField, class GeoMesh >
341 static List<Eigen::MatrixXd> PtrList2EigenBC(
342 PtrList<GeometricField<tensor, PatchField, GeoMesh >> & fields,
343 label Nfields = -1);
344
345 //----------------------------------------------------------------------
357 template<template<class> class PatchField, class GeoMesh >
358 static GeometricField<scalar, PatchField, GeoMesh> Eigen2field(
359 GeometricField<scalar, PatchField, GeoMesh>& field,
360 Eigen::VectorXd& eigen_vector, bool correctBC = true);
361
362 //----------------------------------------------------------------------
374 template<template<class> class PatchField, class GeoMesh>
375 static GeometricField<vector, PatchField, GeoMesh> Eigen2field(
376 GeometricField<vector, PatchField, GeoMesh>& field,
377 Eigen::VectorXd& eigen_vector, bool correctBC = true);
378
379 //----------------------------------------------------------------------
392 template<template<class> class PatchField, class GeoMesh>
393 static GeometricField<scalar, PatchField, GeoMesh> Eigen2field(
394 GeometricField<scalar, PatchField, GeoMesh>& field,
395 Eigen::VectorXd& eigen_vector, List<Eigen::VectorXd>& eigen_vector_boundary);
396
397 //----------------------------------------------------------------------
410 template<template<class> class PatchField, class GeoMesh>
411 static GeometricField<vector, PatchField, GeoMesh> Eigen2field(
412 GeometricField<vector, PatchField, GeoMesh>& field,
413 Eigen::VectorXd& eigen_vector, List<Eigen::VectorXd>& eigen_vector_boundary);
414
415
416 //----------------------------------------------------------------------
428 template<template<class> class PatchField, class GeoMesh>
429 static GeometricField<tensor, PatchField, GeoMesh> Eigen2field(
430 GeometricField<tensor, PatchField, GeoMesh>& field,
431 Eigen::VectorXd& eigen_vector, bool correctBC = true);
432
433 //----------------------------------------------------------------------
444 template <class Type>
445 static Field<Type> Eigen2field(
446 Field<Type>& field, Eigen::MatrixXd& matrix, bool correctBC = true);
447
448 //----------------------------------------------------------------------
461 template<class Type, template<class> class PatchField, class GeoMesh>
462 static std::tuple<Eigen::MatrixXd, Eigen::VectorXd> projectFvMatrix(
463 fvMatrix<Type>& matrix,
464 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes, label Nmodes);
465
466
467 //----------------------------------------------------------------------
493 template<class Type, template<class> class PatchField, class GeoMesh >
494 static Eigen::VectorXd projectField(
495 GeometricField<Type, PatchField, GeoMesh>& field,
496 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes,
497 label Nmodes);
498
499 //----------------------------------------------------------------------
511 template<class Type, template<class> class PatchField, class GeoMesh >
512 static Eigen::MatrixXd MassMatrix(
513 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes, label Nmodes);
514
515 //----------------------------------------------------------------------
527 template<class Type>
528 static std::tuple<List<Eigen::SparseMatrix<double >>, List<Eigen::VectorXd >>
529 LFvMatrix2LSM(PtrList<fvMatrix<Type >> & MatrixList);
530
531 //--------------------------------------------------------------------------
540 template <class type_matrix>
541 static Eigen::Matrix<type_matrix, Eigen::Dynamic, Eigen::Dynamic>
542 List2EigenMatrix ( List<type_matrix> list );
543
544 //--------------------------------------------------------------------------
553 template <class type_matrix>
554 static List<type_matrix> EigenMatrix2List (
555 Eigen::Matrix<type_matrix, Eigen::Dynamic, Eigen::Dynamic> matrix );
556
557 //----------------------------------------------------------------------
566 template <class type_list>
567 static Eigen::MatrixXd field2Eigen(const List<type_list>& list);
568};
569
570#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:59
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:514
static Eigen::VectorXd field2Eigen(GeometricField< tensor, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
Definition Foam2Eigen.C:40
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:245
static Eigen::Map< Eigen::MatrixXd > field2EigenMap(GeometricField< scalar, PatchField, GeoMesh > &field)
Convert a scalar OpenFOAM field into an Eigen Vector.
Definition Foam2Eigen.C:97
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 void fvMatrix2EigenV(fvMatrix< type_foam_matrix > &foam_matrix, type_B &b)
Convert a ldu OpenFOAM matrix into a source vector b.
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:270
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:646
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::Map< Eigen::MatrixXd > field2EigenMapBC(GeometricField< scalar, PatchField, GeoMesh > &field, int BC_index)
Map a scalar OpenFOAM field boundary into an Eigen Matrix.
Definition Foam2Eigen.C:105
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.
volScalarField & A