initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / parallelProcessing / reconstructPar / reconstructPar.C
blobdb14efe3e5eb2c8a395d6815641c0d43bda0e65d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     reconstructPar
28 Description
29     Reconstructs a mesh and fields of a case that is decomposed for parallel
30     execution of OpenFOAM.
32 \*---------------------------------------------------------------------------*/
34 #include "argList.H"
35 #include "timeSelector.H"
37 #include "fvCFD.H"
38 #include "IOobjectList.H"
39 #include "processorMeshes.H"
40 #include "fvFieldReconstructor.H"
41 #include "pointFieldReconstructor.H"
42 #include "reconstructLagrangian.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 int main(int argc, char *argv[])
48     argList::noParallel();
49     timeSelector::addOptions();
50 #   include "addRegionOption.H"
51 #   include "setRootCase.H"
52 #   include "createTime.H"
54     // determine the processor count directly
55     label nProcs = 0;
56     while (dir(args.path()/(word("processor") + name(nProcs))))
57     {
58         ++nProcs;
59     }
61     if (!nProcs)
62     {
63         FatalErrorIn(args.executable())
64             << "No processor* directories found"
65             << exit(FatalError);
66     }
68     // Create the processor databases
69     PtrList<Time> databases(nProcs);
71     forAll (databases, procI)
72     {
73         databases.set
74         (
75             procI,
76             new Time
77             (
78                 Time::controlDictName,
79                 args.rootPath(),
80                 args.caseName()/fileName(word("processor") + name(procI))
81             )
82         );
83     }
85     // use the times list from the master processor
86     // and select a subset based on the command-line options
87     instantList timeDirs = timeSelector::select
88     (
89         databases[0].times(),
90         args
91     );
93     if (!timeDirs.size())
94     {
95         FatalErrorIn(args.executable())
96             << "No times selected"
97             << exit(FatalError);
98     }
100 #   include "createNamedMesh.H"
101     fileName regionPrefix = "";
102     if (regionName != fvMesh::defaultRegion)
103     {
104         regionPrefix = regionName;
105     }
108     // Set all times (on reconstructed mesh and on processor meshes)
109     runTime.setTime(timeDirs[0], 0);
110     mesh.readUpdate();
112     forAll (databases, procI)
113     {
114         databases[procI].setTime(timeDirs[0], 0);
115     }
117     // Read all meshes and addressing to reconstructed mesh
118     processorMeshes procMeshes(databases, regionName);
120     // check face addressing for meshes that have been decomposed
121     // with a very old foam version
122 #   include "checkFaceAddressingComp.H"
124     // Loop over all times
125     forAll (timeDirs, timeI)
126     {
127         // Set time for global database
128         runTime.setTime(timeDirs[timeI], timeI);
130         Info << "Time = " << runTime.timeName() << endl << endl;
132         // Set time for all databases
133         forAll (databases, procI)
134         {
135             databases[procI].setTime(timeDirs[timeI], timeI);
136         }
138         // Check if any new meshes need to be read.
139         fvMesh::readUpdateState meshStat = mesh.readUpdate();
141         fvMesh::readUpdateState procStat = procMeshes.readUpdate();
143         if (procStat == fvMesh::POINTS_MOVED)
144         {
145             // Reconstruct the points for moving mesh cases and write them out
146             procMeshes.reconstructPoints(mesh);
147         }
148         else if (meshStat != procStat)
149         {
150             WarningIn(args.executable())
151                 << "readUpdate for the reconstructed mesh:" << meshStat << nl
152                 << "readUpdate for the processor meshes  :" << procStat << nl
153                 << "These should be equal or your addressing"
154                 << " might be incorrect."
155                 << " Please check your time directories for any "
156                 << "mesh directories." << endl;
157         }
160         // Get list of objects from processor0 database
161         IOobjectList objects(procMeshes.meshes()[0], databases[0].timeName());
164         // If there are any FV fields, reconstruct them
166         if
167         (
168             objects.lookupClass(volScalarField::typeName).size()
169          || objects.lookupClass(volVectorField::typeName).size()
170          || objects.lookupClass(volSphericalTensorField::typeName).size()
171          || objects.lookupClass(volSymmTensorField::typeName).size()
172          || objects.lookupClass(volTensorField::typeName).size()
173          || objects.lookupClass(surfaceScalarField::typeName).size()
174         )
175         {
176             Info << "Reconstructing FV fields" << nl << endl;
178             fvFieldReconstructor fvReconstructor
179             (
180                 mesh,
181                 procMeshes.meshes(),
182                 procMeshes.faceProcAddressing(),
183                 procMeshes.cellProcAddressing(),
184                 procMeshes.boundaryProcAddressing()
185             );
187             fvReconstructor.reconstructFvVolumeFields<scalar>(objects);
188             fvReconstructor.reconstructFvVolumeFields<vector>(objects);
189             fvReconstructor.reconstructFvVolumeFields<sphericalTensor>(objects);
190             fvReconstructor.reconstructFvVolumeFields<symmTensor>(objects);
191             fvReconstructor.reconstructFvVolumeFields<tensor>(objects);
193             fvReconstructor.reconstructFvSurfaceFields<scalar>(objects);
194         }
195         else
196         {
197             Info << "No FV fields" << nl << endl;
198         }
201         // If there are any point fields, reconstruct them
202         if
203         (
204             objects.lookupClass(pointScalarField::typeName).size()
205          || objects.lookupClass(pointVectorField::typeName).size()
206          || objects.lookupClass(pointSphericalTensorField::typeName).size()
207          || objects.lookupClass(pointSymmTensorField::typeName).size()
208          || objects.lookupClass(pointTensorField::typeName).size()
209         )
210         {
211             Info << "Reconstructing point fields" << nl << endl;
213             pointMesh pMesh(mesh);
214             PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
216             forAll (pMeshes, procI)
217             {
218                 pMeshes.set(procI, new pointMesh(procMeshes.meshes()[procI]));
219             }
221             pointFieldReconstructor pointReconstructor
222             (
223                 pMesh,
224                 pMeshes,
225                 procMeshes.pointProcAddressing(),
226                 procMeshes.boundaryProcAddressing()
227             );
229             pointReconstructor.reconstructFields<scalar>(objects);
230             pointReconstructor.reconstructFields<vector>(objects);
231             pointReconstructor.reconstructFields<sphericalTensor>(objects);
232             pointReconstructor.reconstructFields<symmTensor>(objects);
233             pointReconstructor.reconstructFields<tensor>(objects);
234         }
235         else
236         {
237             Info << "No point fields" << nl << endl;
238         }
241         // If there are any clouds, reconstruct them.
242         // The problem is that a cloud of size zero will not get written so
243         // in pass 1 we determine the cloud names and per cloud name the
244         // fields. Note that the fields are stored as IOobjectList from
245         // the first processor that has them. They are in pass2 only used
246         // for name and type (scalar, vector etc).
248         HashTable<IOobjectList> cloudObjects;
250         forAll (databases, procI)
251         {
252             fileNameList cloudDirs
253             (
254                 readDir
255                 (
256                     databases[procI].timePath()/regionPrefix/"lagrangian",
257                     fileName::DIRECTORY
258                 )
259             );
261             forAll (cloudDirs, i)
262             {
263                 // Check if we already have cloud objects for this cloudname.
264                 HashTable<IOobjectList>::const_iterator iter =
265                     cloudObjects.find(cloudDirs[i]);
267                 if (iter == cloudObjects.end())
268                 {
269                     // Do local scan for valid cloud objects.
270                     IOobjectList sprayObjs
271                     (
272                         procMeshes.meshes()[procI],
273                         databases[procI].timeName(),
274                         "lagrangian"/cloudDirs[i]
275                     );
277                     IOobject* positionsPtr = sprayObjs.lookup("positions");
279                     if (positionsPtr)
280                     {
281                         cloudObjects.insert(cloudDirs[i], sprayObjs);
282                     }
283                 }
284             }
285         }
288         if (cloudObjects.size() > 0)
289         {
290             // Pass2: reconstruct the cloud
291             forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
292             {
293                 const word cloudName = string::validate<word>(iter.key());
295                 // Objects (on arbitrary processor)
296                 const IOobjectList& sprayObjs = iter();
298                 Info<< "Reconstructing lagrangian fields for cloud "
299                     << cloudName << nl << endl;
301                 reconstructLagrangianPositions
302                 (
303                     mesh,
304                     cloudName,
305                     procMeshes.meshes(),
306                     procMeshes.faceProcAddressing(),
307                     procMeshes.cellProcAddressing()
308                 );
309                 reconstructLagrangianFields<label>
310                 (
311                     cloudName,
312                     mesh,
313                     procMeshes.meshes(),
314                     sprayObjs
315                 );
316                 reconstructLagrangianFields<scalar>
317                 (
318                     cloudName,
319                     mesh,
320                     procMeshes.meshes(),
321                     sprayObjs
322                 );
323                 reconstructLagrangianFields<vector>
324                 (
325                     cloudName,
326                     mesh,
327                     procMeshes.meshes(),
328                     sprayObjs
329                 );
330                 reconstructLagrangianFields<sphericalTensor>
331                 (
332                     cloudName,
333                     mesh,
334                     procMeshes.meshes(),
335                     sprayObjs
336                 );
337                 reconstructLagrangianFields<symmTensor>
338                 (
339                     cloudName,
340                     mesh,
341                     procMeshes.meshes(),
342                     sprayObjs
343                 );
344                 reconstructLagrangianFields<tensor>
345                 (
346                     cloudName,
347                     mesh,
348                     procMeshes.meshes(),
349                     sprayObjs
350                 );
351             }
352         }
353         else
354         {
355             Info << "No lagrangian fields" << nl << endl;
356         }
358         // If there are any "uniform" directories copy them from
359         // the master processor.
361         fileName uniformDir0 = databases[0].timePath()/"uniform";
362         if (dir(uniformDir0))
363         {
364             cp(uniformDir0, runTime.timePath());
365         }
366     }
368     Info<< "End.\n" << endl;
370     return 0;
374 // ************************************************************************* //