Loading...
Searching...
No Matches
perform_POD.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2Copyright (C) 2017 by the ITHACA-FV authors
3
4License
5 This file is part of ITHACA-FV
6
7 ITHACA-FV is free software: you can redistribute it and/or modify
8 it under the terms of the GNU Lesser General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 ITHACA-FV is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU Lesser General Public License for more details.
16
17 You should have received a copy of the GNU Lesser General Public License
18 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
19
20Class
21 POD
22
23Description
24 Application to perform a POD decomposition of a general field on a given case
25
26SourceFiles
27 perform_POD.C
28
29\*---------------------------------------------------------------------------*/
30
35
38
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "IFstream.H"
43#include "primitiveFields.H"
44#include "FieldFields.H"
45#include "scalarMatrices.H"
46#include "SortableList.H"
47#include "volFieldsFwd.H"
48#include "forces.H"
49#include "forceCoeffs.H"
50#include "volFields.H"
51#include <iostream>
52#include <fstream>
53#include <sstream>
54#include <vector>
55#include <string>
56#include <stdio.h>
57#include "ITHACAPOD.H"
58#include "ITHACAparameters.H"
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61int main(int argc, char *argv[])
62{
63
64#include "setRootCase.H"
65#include "createTime.H"
66#include "createMesh.H"
67
68 PtrList<volVectorField> Vfield;
69 PtrList<volScalarField> Sfield;
70 PtrList<volVectorField> Vmodes;
71 PtrList<volScalarField> Smodes;
72
74 runTime);
75
76 bool pod_exist;
77 struct stat sb;
78
79 if (stat("./ITHACAoutput/POD", &sb) == 0 && S_ISDIR(sb.st_mode))
80 {
81 pod_exist = true;
82 }
83 else
84 {
85 pod_exist = false;
86 Info << "POD don't exist, performing a POD decomposition" << endl;
87 mkDir("./ITHACAoutput/POD");
88 system("ln -s ../../constant ./ITHACAoutput/POD/constant");
89 system("ln -s ../../0 ./ITHACAoutput/POD/0");
90 system("ln -s ../../system ./ITHACAoutput/POD/system");
91 }
92 if(pod_exist == 1)
93 {
94 Info << "The POD has already been performed, please delete the ITHACAoutput folder and try again." << endl;
95 exit(0);
96 }
97
98 // Initialize Variables
99 label nSnapshots;
100 label startTime;
101 label endTime;
102
103 //Read FORCESdict
104 IOdictionary ITHACAPODdict
105 (
106 IOobject
107 (
108 "ITHACAPODdict",
109 runTime.system(),
110 mesh,
111 IOobject::MUST_READ,
112 IOobject::NO_WRITE
113 )
114 );
115
116 // Get times list from the case folder
117 instantList Times = runTime.times();
118
119 // Read Initial and last time from the POD dictionary
120 const entry* existnsnap = ITHACAPODdict.lookupEntryPtr("Nsnapshots", false, true);
121 const entry* existLT = ITHACAPODdict.lookupEntryPtr("FinalTime", false, true);
122
123 // Initiate variable from PODSolverDict
124 if ((existnsnap) && (existLT))
125 {
126 Info << "Error you cannot define LatestTime and NSnapShots together" << endl;
127 abort();
128 }
129 else if (existnsnap)
130 {
131 scalar InitialTime = ITHACAPODdict.lookupOrDefault<scalar>("InitialTime", 0);
132 nSnapshots = readScalar(ITHACAPODdict.lookup("Nsnapshots"));
133 startTime = Time::findClosestTimeIndex(runTime.times(), InitialTime);
134 nSnapshots = min(nSnapshots , Times.size() - startTime);
135 endTime = startTime + nSnapshots - 1;
136 Info << nSnapshots << endl;
137 }
138 else
139 {
140 scalar InitialTime = ITHACAPODdict.lookupOrDefault<scalar>("InitialTime", 0);
141 scalar FinalTime = ITHACAPODdict.lookupOrDefault<scalar>("FinalTime", 100000000000000);
142 endTime = Time::findClosestTimeIndex(runTime.times(), FinalTime);
143 startTime = Time::findClosestTimeIndex(runTime.times(), InitialTime);
144 nSnapshots = endTime - startTime + 1;
145 if (InitialTime > FinalTime)
146 {
147 Info << "FinalTime cannot be smaller than the InitialTime check your ITHACAPODdict file\n" << endl;
148 abort();
149 }
150 }
151 // Print out some Infos
152 Info << "startTime: " << startTime << "\n" << "endTime: " << endTime << "\n" << "nSnapshots: " << nSnapshots << "\n" << endl;
153
154 // Set the initial time
155 runTime.setTime(Times[startTime], startTime);
156
157 wordList fieldlist
158 (
159 ITHACAPODdict.lookup("fields")
160 );
161
162 //word Name = ITHACAPODdict.lookup("fieldName");
163 //word type = ITHACAPODdict.lookup("type");
164
165 if (startTime == endTime)
166 {
167 Info << "The case has no snapshots to process, exiting the code" << endl;
168 exit(0);
169 }
170
171 for (label k = 0; k < fieldlist.size(); k++)
172 {
173 dictionary& subDict = ITHACAPODdict.subDict(fieldlist[k]);
174 scalar nmodes = readScalar(subDict.lookup("nmodes"));
175 word field_name(subDict.lookup("field_name"));
176 word field_type(subDict.lookup("field_type"));
177
178 for (label i = startTime; i < endTime + 1; i++)
179 {
180 Info << "Reading snapshot " << i << " for field " << field_name << endl;
181 runTime.setTime(Times[i], i);
182 mesh.readUpdate();
183
184 if (field_type == "vector")
185 {
186
187 volVectorField vector_field
188 (
189 IOobject
190 (
191 field_name,
192 runTime.timeName(),
193 mesh,
194 IOobject::MUST_READ
195 ),
196 mesh
197 );
198 Vfield.append(vector_field.clone());
199 }
200
201 if (field_type == "scalar")
202 {
203 volScalarField scalar_field
204 (
205 IOobject
206 (
207 field_name,
208 runTime.timeName(),
209 mesh,
210 IOobject::MUST_READ
211 ),
212 mesh
213 );
214 Sfield.append(scalar_field.clone());
215 }
216 }
217
218 if (field_type == "vector")
219 {
220 ITHACAPOD::getModes(Vfield, Vmodes, field_name, 0, 0, 0, nmodes, para->correctBC);
221 }
222 if (field_type == "scalar")
223 {
224 ITHACAPOD::getModes(Sfield, Smodes, field_name, 0, 0, 0, nmodes, para->correctBC);
225 }
226
227 Vfield.clear();
228 Sfield.clear();
229 }
230 Info << endl;
231 Info << "End\n" << endl;
232 return 0;
233}
234
235
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
Definition ITHACAPOD.C:93
int main(int argc, char *argv[])
Definition perform_POD.C:61
label i
Definition pEqn.H:46