Loading...
Searching...
No Matches
ITHACADMD.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 ITHACADMD
26Description
27 Implementation of a DMD of a general field
28SourceFiles
29 ITHACADMD.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ITHACADMD_H
38#define ITHACADMD_H
39
40#include "fvCFD.H"
41#include "ITHACAutilities.H"
42#include "ITHACAstream.H"
43#include "ITHACAparameters.H"
44#include "Foam2Eigen.H"
45#include "EigenFunctions.H"
46#include "ITHACAPOD.H"
47#include <functional>
48#include "Modes.H"
49#pragma GCC diagnostic push
50#pragma GCC diagnostic ignored "-Wold-style-cast"
51#pragma GCC diagnostic ignored "-Wnon-virtual-dtor"
52#include <Spectra/GenEigsSolver.h>
53#include <Spectra/SymEigsSolver.h>
54#include <Eigen/Eigen>
55#include <unsupported/Eigen/SparseExtra>
56#include <unsupported/Eigen/MatrixFunctions>
57#include <redsvd>
58#pragma GCC diagnostic pop
59
60/*---------------------------------------------------------------------------*\
61 Class reductionProblem Declaration
62\*---------------------------------------------------------------------------*/
63
68template<class Type, template<class> class PatchField, class GeoMesh>
70{
71 public:
72
79 ITHACADMD(PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
80 double dt);
81
83 PtrList<GeometricField<Type, PatchField, GeoMesh >> snapshotsDMD;
84
87
90
92 Eigen::VectorXcd eigenValues;
93
95 Eigen::MatrixXcd DMDEigenModes;
96
98 List<Eigen::MatrixXcd> DMDEigenModesBC;
99
101 Eigen::MatrixXcd PODm;
102
104 List<Eigen::MatrixXcd> PODmBC;
105
107 Eigen::VectorXd Amplitudes;
108
110 Eigen::MatrixXcd dynamics;
111
113 label NSnaps;
114
117
120
122 bool redSVD;
123
124 //--------------------------------------------------------------------------
131 void getModes(label SVD_rank = -1, bool exact = true,
132 bool exportDMDmodes = false);
133
134 //--------------------------------------------------------------------------
138
139 //--------------------------------------------------------------------------
144 void exportEigs(word exportFolder);
145
146 //--------------------------------------------------------------------------
153 void getDynamics(double tStart, double tFinal, double dt);
154
155 //--------------------------------------------------------------------------
158 void reconstruct(word exportFolder, word fieldName);
159};
160
163
164#endif
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
ITHACADMD< scalar, fvPatchField, volMesh > ITHACADMDvolScalar
Definition ITHACADMD.H:161
ITHACADMD< vector, fvPatchField, volMesh > ITHACADMDvolVector
Definition ITHACADMD.H:162
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
Header file of the Modes class.
Class of the computation of the DMD, it exploits the SVD methods.
Definition ITHACADMD.H:70
PtrList< GeometricField< scalar, fvPatchField, volMesh > > snapshotsDMD
Definition ITHACADMD.H:83
void convert2Foam()
Convert the EigenModes in Matrix form into OpenFOAM GeometricFields.
Definition ITHACADMD.C:216
ITHACADMD(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, double dt)
Constructs the object.
Definition ITHACADMD.C:37
void getDynamics(double tStart, double tFinal, double dt)
Export the dynamics of DMD given an initial time step, a final one and a time step.
Definition ITHACADMD.C:162
void reconstruct(word exportFolder, word fieldName)
Reconstruct and export the solution using the computed dynamics.
Definition ITHACADMD.C:191
void exportEigs(word exportFolder)
Export the eigenvalues in numpy format.
Definition ITHACADMD.C:180
Modes< scalar, fvPatchField, volMesh > DMDmodesImag
Definition ITHACADMD.H:89
List< Eigen::MatrixXcd > DMDEigenModesBC
Definition ITHACADMD.H:98
void getModes(label SVD_rank=-1, bool exact=true, bool exportDMDmodes=false)
Get the DMD modes.
Definition ITHACADMD.C:50
Modes< scalar, fvPatchField, volMesh > DMDmodesReal
Definition ITHACADMD.H:86
Implementation of a container class derived from PtrList.
Definition Modes.H:69