1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 Reconstructs a mesh and fields of a case that is decomposed for parallel
29 execution of OpenFOAM.
31 \*---------------------------------------------------------------------------*/
34 #include "timeSelector.H"
37 #include "IOobjectList.H"
38 #include "processorMeshes.H"
39 #include "fvFieldReconstructor.H"
40 #include "pointFieldReconstructor.H"
41 #include "reconstructLagrangian.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 int main(int argc, char *argv[])
47 // enable -constant ... if someone really wants it
48 // enable -zeroTime to prevent accidentally trashing the initial fields
49 timeSelector::addOptions(true, true);
50 argList::noParallel();
51 # include "addRegionOption.H"
56 "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
57 "regular expressions not currently supported"
63 "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
64 "regular expressions not currently supported, "
65 "positions always included."
67 argList::addBoolOption
70 "skip reconstructing lagrangian positions and fields"
73 # include "setRootCase.H"
74 # include "createTime.H"
76 HashSet<word> selectedFields;
77 if (args.optionFound("fields"))
79 args.optionLookup("fields")() >> selectedFields;
82 const bool noLagrangian = args.optionFound("noLagrangian");
84 HashSet<word> selectedLagrangianFields;
85 if (args.optionFound("lagrangianFields"))
89 FatalErrorIn(args.executable())
90 << "Cannot specify noLagrangian and lagrangianFields "
91 << "options together."
95 args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
98 // determine the processor count directly
100 while (isDir(args.path()/(word("processor") + name(nProcs))))
107 FatalErrorIn(args.executable())
108 << "No processor* directories found"
112 // Create the processor databases
113 PtrList<Time> databases(nProcs);
115 forAll(databases, procI)
122 Time::controlDictName,
124 args.caseName()/fileName(word("processor") + name(procI))
129 // use the times list from the master processor
130 // and select a subset based on the command-line options
131 instantList timeDirs = timeSelector::select
133 databases[0].times(),
137 if (timeDirs.empty())
139 FatalErrorIn(args.executable())
140 << "No times selected"
144 # include "createNamedMesh.H"
145 word regionDir = word::null;
146 if (regionName != fvMesh::defaultRegion)
148 regionDir = regionName;
151 // Set all times on processor meshes equal to reconstructed mesh
152 forAll(databases, procI)
154 databases[procI].setTime(runTime.timeName(), runTime.timeIndex());
157 // Read all meshes and addressing to reconstructed mesh
158 processorMeshes procMeshes(databases, regionName);
161 // check face addressing for meshes that have been decomposed
162 // with a very old foam version
163 # include "checkFaceAddressingComp.H"
165 // Loop over all times
166 forAll(timeDirs, timeI)
168 // Set time for global database
169 runTime.setTime(timeDirs[timeI], timeI);
171 Info<< "Time = " << runTime.timeName() << endl << endl;
173 // Set time for all databases
174 forAll(databases, procI)
176 databases[procI].setTime(timeDirs[timeI], timeI);
179 // Check if any new meshes need to be read.
180 fvMesh::readUpdateState meshStat = mesh.readUpdate();
182 fvMesh::readUpdateState procStat = procMeshes.readUpdate();
184 if (procStat == fvMesh::POINTS_MOVED)
186 // Reconstruct the points for moving mesh cases and write them out
187 procMeshes.reconstructPoints(mesh);
189 else if (meshStat != procStat)
191 WarningIn(args.executable())
192 << "readUpdate for the reconstructed mesh:" << meshStat << nl
193 << "readUpdate for the processor meshes :" << procStat << nl
194 << "These should be equal or your addressing"
195 << " might be incorrect."
196 << " Please check your time directories for any "
197 << "mesh directories." << endl;
201 // Get list of objects from processor0 database
202 IOobjectList objects(procMeshes.meshes()[0], databases[0].timeName());
205 // If there are any FV fields, reconstruct them
206 Info<< "Reconstructing FV fields" << nl << endl;
208 fvFieldReconstructor fvReconstructor
212 procMeshes.faceProcAddressing(),
213 procMeshes.cellProcAddressing(),
214 procMeshes.boundaryProcAddressing()
217 fvReconstructor.reconstructFvVolumeInternalFields<scalar>
222 fvReconstructor.reconstructFvVolumeInternalFields<vector>
227 fvReconstructor.reconstructFvVolumeInternalFields<sphericalTensor>
232 fvReconstructor.reconstructFvVolumeInternalFields<symmTensor>
237 fvReconstructor.reconstructFvVolumeInternalFields<tensor>
243 fvReconstructor.reconstructFvVolumeFields<scalar>
248 fvReconstructor.reconstructFvVolumeFields<vector>
253 fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
258 fvReconstructor.reconstructFvVolumeFields<symmTensor>
263 fvReconstructor.reconstructFvVolumeFields<tensor>
269 fvReconstructor.reconstructFvSurfaceFields<scalar>
274 fvReconstructor.reconstructFvSurfaceFields<vector>
279 fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
284 fvReconstructor.reconstructFvSurfaceFields<symmTensor>
289 fvReconstructor.reconstructFvSurfaceFields<tensor>
295 if (fvReconstructor.nReconstructed() == 0)
297 Info<< "No FV fields" << nl << endl;
302 Info<< "Reconstructing point fields" << nl << endl;
304 pointMesh pMesh(mesh);
305 PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
307 forAll(pMeshes, procI)
309 pMeshes.set(procI, new pointMesh(procMeshes.meshes()[procI]));
312 pointFieldReconstructor pointReconstructor
316 procMeshes.pointProcAddressing(),
317 procMeshes.boundaryProcAddressing()
320 pointReconstructor.reconstructFields<scalar>
325 pointReconstructor.reconstructFields<vector>
330 pointReconstructor.reconstructFields<sphericalTensor>
335 pointReconstructor.reconstructFields<symmTensor>
340 pointReconstructor.reconstructFields<tensor>
346 if (pointReconstructor.nReconstructed() == 0)
348 Info<< "No point fields" << nl << endl;
353 // If there are any clouds, reconstruct them.
354 // The problem is that a cloud of size zero will not get written so
355 // in pass 1 we determine the cloud names and per cloud name the
356 // fields. Note that the fields are stored as IOobjectList from
357 // the first processor that has them. They are in pass2 only used
358 // for name and type (scalar, vector etc).
362 HashTable<IOobjectList> cloudObjects;
364 forAll(databases, procI)
366 fileNameList cloudDirs
370 databases[procI].timePath() / regionDir / cloud::prefix,
377 // Check if we already have cloud objects for this cloudname
378 HashTable<IOobjectList>::const_iterator iter =
379 cloudObjects.find(cloudDirs[i]);
381 if (iter == cloudObjects.end())
383 // Do local scan for valid cloud objects
384 IOobjectList sprayObjs
386 procMeshes.meshes()[procI],
387 databases[procI].timeName(),
388 cloud::prefix/cloudDirs[i]
391 IOobject* positionsPtr = sprayObjs.lookup("positions");
395 cloudObjects.insert(cloudDirs[i], sprayObjs);
402 if (cloudObjects.size())
404 // Pass2: reconstruct the cloud
405 forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
407 const word cloudName = string::validate<word>(iter.key());
409 // Objects (on arbitrary processor)
410 const IOobjectList& sprayObjs = iter();
412 Info<< "Reconstructing lagrangian fields for cloud "
413 << cloudName << nl << endl;
415 reconstructLagrangianPositions
420 procMeshes.faceProcAddressing(),
421 procMeshes.cellProcAddressing()
423 reconstructLagrangianFields<label>
429 selectedLagrangianFields
431 reconstructLagrangianFieldFields<label>
437 selectedLagrangianFields
439 reconstructLagrangianFields<scalar>
445 selectedLagrangianFields
447 reconstructLagrangianFieldFields<scalar>
453 selectedLagrangianFields
455 reconstructLagrangianFields<vector>
461 selectedLagrangianFields
463 reconstructLagrangianFieldFields<vector>
469 selectedLagrangianFields
471 reconstructLagrangianFields<sphericalTensor>
477 selectedLagrangianFields
479 reconstructLagrangianFieldFields<sphericalTensor>
485 selectedLagrangianFields
487 reconstructLagrangianFields<symmTensor>
493 selectedLagrangianFields
495 reconstructLagrangianFieldFields<symmTensor>
501 selectedLagrangianFields
503 reconstructLagrangianFields<tensor>
509 selectedLagrangianFields
511 reconstructLagrangianFieldFields<tensor>
517 selectedLagrangianFields
523 Info<< "No lagrangian fields" << nl << endl;
527 // If there are any "uniform" directories copy them from
528 // the master processor
530 fileName uniformDir0 = databases[0].timePath()/"uniform";
531 if (isDir(uniformDir0))
533 cp(uniformDir0, runTime.timePath());
537 Info<< "End.\n" << endl;
543 // ************************************************************************* //