BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / parallelProcessing / reconstructPar / reconstructPar.C
blob37d4a04f2221e72a6f59a7c505062bf4f3678784
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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/>.
24 Application
25     reconstructPar
27 Description
28     Reconstructs a mesh and fields of a case that is decomposed for parallel
29     execution of OpenFOAM.
31 \*---------------------------------------------------------------------------*/
33 #include "argList.H"
34 #include "timeSelector.H"
36 #include "fvCFD.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"
52     argList::addOption
53     (
54         "fields",
55         "list",
56         "specify a list of fields to be reconstructed. Eg, '(U T p)' - "
57         "regular expressions not currently supported"
58     );
59     argList::addOption
60     (
61         "lagrangianFields",
62         "list",
63         "specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
64         "regular expressions not currently supported, "
65         "positions always included."
66     );
67     argList::addBoolOption
68     (
69         "noLagrangian",
70         "skip reconstructing lagrangian positions and fields"
71     );
73 #   include "setRootCase.H"
74 #   include "createTime.H"
76     HashSet<word> selectedFields;
77     if (args.optionFound("fields"))
78     {
79         args.optionLookup("fields")() >> selectedFields;
80     }
82     const bool noLagrangian = args.optionFound("noLagrangian");
84     HashSet<word> selectedLagrangianFields;
85     if (args.optionFound("lagrangianFields"))
86     {
87         if (noLagrangian)
88         {
89             FatalErrorIn(args.executable())
90                 << "Cannot specify noLagrangian and lagrangianFields "
91                 << "options together."
92                 << exit(FatalError);
93         }
95         args.optionLookup("lagrangianFields")() >> selectedLagrangianFields;
96     }
98     // determine the processor count directly
99     label nProcs = 0;
100     while (isDir(args.path()/(word("processor") + name(nProcs))))
101     {
102         ++nProcs;
103     }
105     if (!nProcs)
106     {
107         FatalErrorIn(args.executable())
108             << "No processor* directories found"
109             << exit(FatalError);
110     }
112     // Create the processor databases
113     PtrList<Time> databases(nProcs);
115     forAll(databases, procI)
116     {
117         databases.set
118         (
119             procI,
120             new Time
121             (
122                 Time::controlDictName,
123                 args.rootPath(),
124                 args.caseName()/fileName(word("processor") + name(procI))
125             )
126         );
127     }
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
132     (
133         databases[0].times(),
134         args
135     );
137     if (timeDirs.empty())
138     {
139         FatalErrorIn(args.executable())
140             << "No times selected"
141             << exit(FatalError);
142     }
144 #   include "createNamedMesh.H"
145     word regionDir = word::null;
146     if (regionName != fvMesh::defaultRegion)
147     {
148         regionDir = regionName;
149     }
151     // Set all times on processor meshes equal to reconstructed mesh
152     forAll(databases, procI)
153     {
154         databases[procI].setTime(runTime.timeName(), runTime.timeIndex());
155     }
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)
167     {
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)
175         {
176             databases[procI].setTime(timeDirs[timeI], timeI);
177         }
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)
185         {
186             // Reconstruct the points for moving mesh cases and write them out
187             procMeshes.reconstructPoints(mesh);
188         }
189         else if (meshStat != procStat)
190         {
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;
198         }
201         // Get list of objects from processor0 database
202         IOobjectList objects(procMeshes.meshes()[0], databases[0].timeName());
204         {
205             // If there are any FV fields, reconstruct them
206             Info<< "Reconstructing FV fields" << nl << endl;
208             fvFieldReconstructor fvReconstructor
209             (
210                 mesh,
211                 procMeshes.meshes(),
212                 procMeshes.faceProcAddressing(),
213                 procMeshes.cellProcAddressing(),
214                 procMeshes.boundaryProcAddressing()
215             );
217             fvReconstructor.reconstructFvVolumeInternalFields<scalar>
218             (
219                 objects,
220                 selectedFields
221             );
222             fvReconstructor.reconstructFvVolumeInternalFields<vector>
223             (
224                 objects,
225                 selectedFields
226             );
227             fvReconstructor.reconstructFvVolumeInternalFields<sphericalTensor>
228             (
229                 objects,
230                 selectedFields
231             );
232             fvReconstructor.reconstructFvVolumeInternalFields<symmTensor>
233             (
234                 objects,
235                 selectedFields
236             );
237             fvReconstructor.reconstructFvVolumeInternalFields<tensor>
238             (
239                 objects,
240                 selectedFields
241             );
243             fvReconstructor.reconstructFvVolumeFields<scalar>
244             (
245                 objects,
246                 selectedFields
247             );
248             fvReconstructor.reconstructFvVolumeFields<vector>
249             (
250                 objects,
251                 selectedFields
252             );
253             fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
254             (
255                 objects,
256                 selectedFields
257             );
258             fvReconstructor.reconstructFvVolumeFields<symmTensor>
259             (
260                 objects,
261                 selectedFields
262             );
263             fvReconstructor.reconstructFvVolumeFields<tensor>
264             (
265                 objects,
266                 selectedFields
267             );
269             fvReconstructor.reconstructFvSurfaceFields<scalar>
270             (
271                 objects,
272                 selectedFields
273             );
274             fvReconstructor.reconstructFvSurfaceFields<vector>
275             (
276                 objects,
277                 selectedFields
278             );
279             fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
280             (
281                 objects,
282                 selectedFields
283             );
284             fvReconstructor.reconstructFvSurfaceFields<symmTensor>
285             (
286                 objects,
287                 selectedFields
288             );
289             fvReconstructor.reconstructFvSurfaceFields<tensor>
290             (
291                 objects,
292                 selectedFields
293             );
295             if (fvReconstructor.nReconstructed() == 0)
296             {
297                 Info<< "No FV fields" << nl << endl;
298             }
299         }
301         {
302             Info<< "Reconstructing point fields" << nl << endl;
304             pointMesh pMesh(mesh);
305             PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
307             forAll(pMeshes, procI)
308             {
309                 pMeshes.set(procI, new pointMesh(procMeshes.meshes()[procI]));
310             }
312             pointFieldReconstructor pointReconstructor
313             (
314                 pMesh,
315                 pMeshes,
316                 procMeshes.pointProcAddressing(),
317                 procMeshes.boundaryProcAddressing()
318             );
320             pointReconstructor.reconstructFields<scalar>
321             (
322                 objects,
323                 selectedFields
324             );
325             pointReconstructor.reconstructFields<vector>
326             (
327                 objects,
328                 selectedFields
329             );
330             pointReconstructor.reconstructFields<sphericalTensor>
331             (
332                 objects,
333                 selectedFields
334             );
335             pointReconstructor.reconstructFields<symmTensor>
336             (
337                 objects,
338                 selectedFields
339             );
340             pointReconstructor.reconstructFields<tensor>
341             (
342                 objects,
343                 selectedFields
344             );
346             if (pointReconstructor.nReconstructed() == 0)
347             {
348                 Info<< "No point fields" << nl << endl;
349             }
350         }
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).
360         if (!noLagrangian)
361         {
362             HashTable<IOobjectList> cloudObjects;
364             forAll(databases, procI)
365             {
366                 fileNameList cloudDirs
367                 (
368                     readDir
369                     (
370                         databases[procI].timePath() / regionDir / cloud::prefix,
371                         fileName::DIRECTORY
372                     )
373                 );
375                 forAll(cloudDirs, i)
376                 {
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())
382                     {
383                         // Do local scan for valid cloud objects
384                         IOobjectList sprayObjs
385                         (
386                             procMeshes.meshes()[procI],
387                             databases[procI].timeName(),
388                             cloud::prefix/cloudDirs[i]
389                         );
391                         IOobject* positionsPtr = sprayObjs.lookup("positions");
393                         if (positionsPtr)
394                         {
395                             cloudObjects.insert(cloudDirs[i], sprayObjs);
396                         }
397                     }
398                 }
399             }
402             if (cloudObjects.size())
403             {
404                 // Pass2: reconstruct the cloud
405                 forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
406                 {
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
416                     (
417                         mesh,
418                         cloudName,
419                         procMeshes.meshes(),
420                         procMeshes.faceProcAddressing(),
421                         procMeshes.cellProcAddressing()
422                     );
423                     reconstructLagrangianFields<label>
424                     (
425                         cloudName,
426                         mesh,
427                         procMeshes.meshes(),
428                         sprayObjs,
429                         selectedLagrangianFields
430                     );
431                     reconstructLagrangianFieldFields<label>
432                     (
433                         cloudName,
434                         mesh,
435                         procMeshes.meshes(),
436                         sprayObjs,
437                         selectedLagrangianFields
438                     );
439                     reconstructLagrangianFields<scalar>
440                     (
441                         cloudName,
442                         mesh,
443                         procMeshes.meshes(),
444                         sprayObjs,
445                         selectedLagrangianFields
446                     );
447                     reconstructLagrangianFieldFields<scalar>
448                     (
449                         cloudName,
450                         mesh,
451                         procMeshes.meshes(),
452                         sprayObjs,
453                         selectedLagrangianFields
454                     );
455                     reconstructLagrangianFields<vector>
456                     (
457                         cloudName,
458                         mesh,
459                         procMeshes.meshes(),
460                         sprayObjs,
461                         selectedLagrangianFields
462                     );
463                     reconstructLagrangianFieldFields<vector>
464                     (
465                         cloudName,
466                         mesh,
467                         procMeshes.meshes(),
468                         sprayObjs,
469                         selectedLagrangianFields
470                     );
471                     reconstructLagrangianFields<sphericalTensor>
472                     (
473                         cloudName,
474                         mesh,
475                         procMeshes.meshes(),
476                         sprayObjs,
477                         selectedLagrangianFields
478                     );
479                     reconstructLagrangianFieldFields<sphericalTensor>
480                     (
481                         cloudName,
482                         mesh,
483                         procMeshes.meshes(),
484                         sprayObjs,
485                         selectedLagrangianFields
486                     );
487                     reconstructLagrangianFields<symmTensor>
488                     (
489                         cloudName,
490                         mesh,
491                         procMeshes.meshes(),
492                         sprayObjs,
493                         selectedLagrangianFields
494                     );
495                     reconstructLagrangianFieldFields<symmTensor>
496                     (
497                         cloudName,
498                         mesh,
499                         procMeshes.meshes(),
500                         sprayObjs,
501                         selectedLagrangianFields
502                     );
503                     reconstructLagrangianFields<tensor>
504                     (
505                         cloudName,
506                         mesh,
507                         procMeshes.meshes(),
508                         sprayObjs,
509                         selectedLagrangianFields
510                     );
511                     reconstructLagrangianFieldFields<tensor>
512                     (
513                         cloudName,
514                         mesh,
515                         procMeshes.meshes(),
516                         sprayObjs,
517                         selectedLagrangianFields
518                     );
519                 }
520             }
521             else
522             {
523                 Info<< "No lagrangian fields" << nl << endl;
524             }
525         }
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))
532         {
533             cp(uniformDir0, runTime.timePath());
534         }
535     }
537     Info<< "End.\n" << endl;
539     return 0;
543 // ************************************************************************* //