BUGFIX: Illegal use of uninitialised value (backport)
[foam-extend-3.2.git] / applications / utilities / parallelProcessing / decomposePar / decomposePar.C
blob192c29bbe97179c693e9507825934d8228473880
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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     decomposePar
28 Description
29     Automatically decomposes a mesh and fields of a case for parallel
30     execution of OpenFOAM.
32 Usage
34     - decomposePar [OPTION]
36     @param -cellDist \n
37     Write the cell distribution as a labelList for use with 'manual'
38     decomposition method and as a volScalarField for post-processing.
40     @param -region regionName \n
41     Decompose named region. Does not check for existence of processor*.
43     @param -copyUniform \n
44     Copy any @a uniform directories too.
46     @param -fields \n
47     Use existing geometry decomposition and convert fields only.
49     @param -filterPatches \n
50     Remove empty patches when decomposing the geometry.
52     @param -force \n
53     Remove any existing @a processor subdirectories before decomposing the
54     geometry.
56     @param -ifRequired \n
57     Only decompose the geometry if the number of domains has changed from a
58     previous decomposition. No @a processor subdirectories will be removed
59     unless the @a -force option is also specified. This option can be used
60     to avoid redundant geometry decomposition (eg, in scripts), but should
61     be used with caution when the underlying (serial) geometry or the
62     decomposition method etc. have been changed between decompositions.
64 \*---------------------------------------------------------------------------*/
66 #include "OSspecific.H"
67 #include "fvCFD.H"
68 #include "IOobjectList.H"
69 #include "processorFvPatchFields.H"
70 #include "domainDecomposition.H"
71 #include "labelIOField.H"
72 #include "scalarIOField.H"
73 #include "vectorIOField.H"
74 #include "sphericalTensorIOField.H"
75 #include "symmTensorIOField.H"
76 #include "tensorIOField.H"
78 #include "tetPointFields.H"
79 #include "elementFields.H"
80 #include "tetFemMatrices.H"
81 #include "tetPointFieldDecomposer.H"
83 #include "pointFields.H"
85 #include "readFields.H"
86 #include "fvFieldDecomposer.H"
87 #include "pointFieldDecomposer.H"
88 #include "lagrangianFieldDecomposer.H"
90 #include "faCFD.H"
91 #include "emptyFaPatch.H"
92 #include "faMeshDecomposition.H"
93 #include "faFieldDecomposer.H"
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 int main(int argc, char *argv[])
99     argList::noParallel();
100 #   include "addRegionOption.H"
101     argList::validOptions.insert("cellDist", "");
102     argList::validOptions.insert("copyUniform", "");
103     argList::validOptions.insert("fields", "");
104     argList::validOptions.insert("filterPatches", "");
105     argList::validOptions.insert("force", "");
106     argList::validOptions.insert("ifRequired", "");
108 #   include "setRootCase.H"
110     word regionName = fvMesh::defaultRegion;
111     word regionDir = word::null;
113     if (args.optionFound("region"))
114     {
115         regionName = args.option("region");
116         regionDir = regionName;
117         Info<< "Decomposing mesh " << regionName << nl << endl;
118     }
121     bool writeCellDist = args.optionFound("cellDist");
122     bool copyUniform = args.optionFound("copyUniform");
123     bool decomposeFieldsOnly = args.optionFound("fields");
124     bool filterPatches = args.optionFound("filterPatches");
125     bool forceOverwrite = args.optionFound("force");
126     bool ifRequiredDecomposition = args.optionFound("ifRequired");
128 #   include "createTime.H"
130     Info<< "Time = " << runTime.timeName() << endl;
132     // Determine the existing processor count directly
133     label nProcs = 0;
134     while
135     (
136         isDir
137         (
138             runTime.path()
139            /(word("processor") + name(nProcs))
140            /runTime.constant()
141            /regionDir
142            /polyMesh::meshSubDir
143         )
144     )
145     {
146         ++nProcs;
147     }
149     // get requested numberOfSubdomains
150     label nDomains = 0;
151     {
152         IOdictionary decompDict
153         (
154             IOobject
155             (
156                 "decomposeParDict",
157                 runTime.time().system(),
158                 regionDir,          // use region if non-standard
159                 runTime,
160                 IOobject::MUST_READ,
161                 IOobject::NO_WRITE,
162                 false
163             )
164         );
166         decompDict.lookup("numberOfSubdomains") >> nDomains;
167     }
169     if (decomposeFieldsOnly)
170     {
171         // Sanity check on previously decomposed case
172         if (nProcs != nDomains)
173         {
174             FatalErrorIn(args.executable())
175                 << "Specified -fields, but the case was decomposed with "
176                 << nProcs << " domains"
177                 << nl
178                 << "instead of " << nDomains
179                 << " domains as specified in decomposeParDict"
180                 << nl
181                 << exit(FatalError);
182         }
183     }
184     else if (nProcs)
185     {
186         bool procDirsProblem = true;
188         if (regionName != fvMesh::defaultRegion)
189         {
190             decomposeFieldsOnly = false;
191             procDirsProblem = false;
192         }
195         if (ifRequiredDecomposition && nProcs == nDomains)
196         {
197             // we can reuse the decomposition
198             decomposeFieldsOnly = true;
199             procDirsProblem = false;
200             forceOverwrite = false;
202             Info<< "Using existing processor directories" << nl;
203         }
205         if (forceOverwrite)
206         {
207             Info<< "Removing " << nProcs
208                 << " existing processor directories" << endl;
210             // remove existing processor dirs
211             // reverse order to avoid gaps if someone interrupts the process
212             for (label procI = nProcs-1; procI >= 0; --procI)
213             {
214                 fileName procDir
215                 (
216                     runTime.path()/(word("processor") + name(procI))
217                 );
219                 rmDir(procDir);
220             }
222             procDirsProblem = false;
223         }
225         if (procDirsProblem)
226         {
227             FatalErrorIn(args.executable())
228                 << "Case is already decomposed with " << nProcs
229                 << " domains, use the -force option or manually" << nl
230                 << "remove processor directories before decomposing. e.g.,"
231                 << nl
232                 << "    rm -rf " << runTime.path().c_str() << "/processor*"
233                 << nl
234                 << exit(FatalError);
235         }
236     }
238     Info<< "Create mesh for region " << regionName << endl;
239     domainDecomposition mesh
240     (
241         IOobject
242         (
243             regionName,
244             runTime.timeName(),
245             runTime
246         )
247     );
249     // Decompose the mesh
250     if (!decomposeFieldsOnly)
251     {
252         mesh.decomposeMesh(filterPatches);
254         mesh.writeDecomposition();
256         if (writeCellDist)
257         {
258             const labelList& procIds = mesh.cellToProc();
260             // Write the decomposition as labelList for use with 'manual'
261             // decomposition method.
262             labelIOList cellDecomposition
263             (
264                 IOobject
265                 (
266                     "cellDecomposition",
267                     mesh.facesInstance(),
268                     mesh,
269                     IOobject::NO_READ,
270                     IOobject::NO_WRITE,
271                     false
272                 ),
273                 procIds
274             );
275             cellDecomposition.write();
277             Info<< nl << "Wrote decomposition to "
278                 << cellDecomposition.objectPath()
279                 << " for use in manual decomposition." << endl;
281             // Write as volScalarField for post-processing
282             volScalarField cellDist
283             (
284                 IOobject
285                 (
286                     "cellDist",
287                     runTime.timeName(),
288                     mesh,
289                     IOobject::NO_READ,
290                     IOobject::NO_WRITE
291                 ),
292                 mesh,
293                 dimensionedScalar("cellDist", dimless, 0),
294                 zeroGradientFvPatchScalarField::typeName
295             );
297             forAll(procIds, celli)
298             {
299                cellDist[celli] = procIds[celli];
300             }
302             cellDist.write();
304             Info<< nl << "Wrote decomposition as volScalarField to "
305                 << cellDist.name() << " for use in post-processing."
306                 << endl;
307         }
308     }
311     // Search for list of objects for this time
312     IOobjectList objects(mesh, runTime.timeName());
314     // Construct the vol fields
315     // ~~~~~~~~~~~~~~~~~~~~~~~~
316     PtrList<volScalarField> volScalarFields;
317     readFields(mesh, objects, volScalarFields);
319     PtrList<volVectorField> volVectorFields;
320     readFields(mesh, objects, volVectorFields);
322     PtrList<volSphericalTensorField> volSphericalTensorFields;
323     readFields(mesh, objects, volSphericalTensorFields);
325     PtrList<volSymmTensorField> volSymmTensorFields;
326     readFields(mesh, objects, volSymmTensorFields);
328     PtrList<volTensorField> volTensorFields;
329     readFields(mesh, objects, volTensorFields);
332     // Construct the surface fields
333     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
334     PtrList<surfaceScalarField> surfaceScalarFields;
335     readFields(mesh, objects, surfaceScalarFields);
336     PtrList<surfaceVectorField> surfaceVectorFields;
337     readFields(mesh, objects, surfaceVectorFields);
338     PtrList<surfaceSphericalTensorField> surfaceSphericalTensorFields;
339     readFields(mesh, objects, surfaceSphericalTensorFields);
340     PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
341     readFields(mesh, objects, surfaceSymmTensorFields);
342     PtrList<surfaceTensorField> surfaceTensorFields;
343     readFields(mesh, objects, surfaceTensorFields);
346     // Construct the point fields
347     // ~~~~~~~~~~~~~~~~~~~~~~~~~~
348     pointMesh pMesh(mesh);
350     PtrList<pointScalarField> pointScalarFields;
351     readFields(pMesh, objects, pointScalarFields);
353     PtrList<pointVectorField> pointVectorFields;
354     readFields(pMesh, objects, pointVectorFields);
356     PtrList<pointSphericalTensorField> pointSphericalTensorFields;
357     readFields(pMesh, objects, pointSphericalTensorFields);
359     PtrList<pointSymmTensorField> pointSymmTensorFields;
360     readFields(pMesh, objects, pointSymmTensorFields);
362     PtrList<pointTensorField> pointTensorFields;
363     readFields(pMesh, objects, pointTensorFields);
366     // Construct the tetPoint fields
367     // ~~~~~~~~~~~~~~~~~~~~~~~~~~
368     tetPolyMesh* tetMeshPtr = NULL;
370     PtrList<tetPointScalarField> tetPointScalarFields;
371     PtrList<tetPointVectorField> tetPointVectorFields;
372     PtrList<tetPointSphericalTensorField> tetPointSphericalTensorFields;
373     PtrList<tetPointSymmTensorField> tetPointSymmTensorFields;
374     PtrList<tetPointTensorField> tetPointTensorFields;
376     PtrList<elementScalarField> elementScalarFields;
377     PtrList<elementVectorField> elementVectorFields;
379     if
380     (
381         objects.lookupClass("tetPointScalarField").size() > 0
382      || objects.lookupClass("tetPointVectorField").size() > 0
383      || objects.lookupClass("tetPointSphericalTensorField").size() > 0
384      || objects.lookupClass("tetPointSymmTensorField").size() > 0
385      || objects.lookupClass("tetPointTensorField").size() > 0
387      || objects.lookupClass("elementScalarField").size() > 0
388      || objects.lookupClass("elementVectorField").size() > 0
389     )
390     {
391         tetMeshPtr = new tetPolyMesh(mesh);
392         tetPolyMesh& tetMesh = *tetMeshPtr;
394         readFields(tetMesh, objects, tetPointScalarFields);
395         readFields(tetMesh, objects, tetPointVectorFields);
396         readFields(tetMesh, objects, tetPointSphericalTensorFields);
397         readFields(tetMesh, objects, tetPointSymmTensorFields);
398         readFields(tetMesh, objects, tetPointTensorFields);
400         readFields(tetMesh, objects, elementScalarFields);
401         readFields(tetMesh, objects, elementVectorFields);
402     }
405     // Construct the Lagrangian fields
406     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
408     fileNameList cloudDirs
409     (
410         readDir(runTime.timePath()/cloud::prefix, fileName::DIRECTORY)
411     );
413     // Particles
414     PtrList<Cloud<indexedParticle> > lagrangianPositions(cloudDirs.size());
415     // Particles per cell
416     PtrList< List<SLList<indexedParticle*>*> > cellParticles(cloudDirs.size());
418     PtrList<PtrList<labelIOField> > lagrangianLabelFields(cloudDirs.size());
419     PtrList<PtrList<scalarIOField> > lagrangianScalarFields(cloudDirs.size());
420     PtrList<PtrList<vectorIOField> > lagrangianVectorFields(cloudDirs.size());
421     PtrList<PtrList<sphericalTensorIOField> > lagrangianSphericalTensorFields
422     (
423         cloudDirs.size()
424     );
425     PtrList<PtrList<symmTensorIOField> > lagrangianSymmTensorFields
426     (
427         cloudDirs.size()
428     );
429     PtrList<PtrList<tensorIOField> > lagrangianTensorFields
430     (
431         cloudDirs.size()
432     );
434     label cloudI = 0;
436     forAll(cloudDirs, i)
437     {
438         IOobjectList sprayObjs
439         (
440             mesh,
441             runTime.timeName(),
442             cloud::prefix/cloudDirs[i]
443         );
445         IOobject* positionsPtr = sprayObjs.lookup("positions");
447         if (positionsPtr)
448         {
449             // Read lagrangian particles
450             // ~~~~~~~~~~~~~~~~~~~~~~~~~
452             Info<< "Identified lagrangian data set: " << cloudDirs[i] << endl;
454             lagrangianPositions.set
455             (
456                 cloudI,
457                 new Cloud<indexedParticle>
458                 (
459                     mesh,
460                     cloudDirs[i],
461                     false
462                 )
463             );
466             // Sort particles per cell
467             // ~~~~~~~~~~~~~~~~~~~~~~~
469             cellParticles.set
470             (
471                 cloudI,
472                 new List<SLList<indexedParticle*>*>
473                 (
474                     mesh.nCells(),
475                     static_cast<SLList<indexedParticle*>*>(NULL)
476                 )
477             );
479             label i = 0;
481             forAllIter
482             (
483                 Cloud<indexedParticle>,
484                 lagrangianPositions[cloudI],
485                 iter
486             )
487             {
488                 iter().index() = i++;
490                 label celli = iter().cell();
492                 // Check
493                 if (celli < 0 || celli >= mesh.nCells())
494                 {
495                     FatalErrorIn(args.executable())
496                         << "Illegal cell number " << celli
497                         << " for particle with index " << iter().index()
498                         << " at position " << iter().position() << nl
499                         << "Cell number should be between 0 and "
500                         << mesh.nCells()-1 << nl
501                         << "On this mesh the particle should be in cell "
502                         << mesh.findCell(iter().position())
503                         << exit(FatalError);
504                 }
506                 if (!cellParticles[cloudI][celli])
507                 {
508                     cellParticles[cloudI][celli] =
509                         new SLList<indexedParticle*>();
510                 }
512                 cellParticles[cloudI][celli]->append(&iter());
513             }
515             // Read fields
516             // ~~~~~~~~~~~
518             IOobjectList lagrangianObjects
519             (
520                 mesh,
521                 runTime.timeName(),
522                 cloud::prefix/cloudDirs[cloudI]
523             );
525             lagrangianFieldDecomposer::readFields
526             (
527                 cloudI,
528                 lagrangianObjects,
529                 lagrangianLabelFields
530             );
532             lagrangianFieldDecomposer::readFields
533             (
534                 cloudI,
535                 lagrangianObjects,
536                 lagrangianScalarFields
537             );
539             lagrangianFieldDecomposer::readFields
540             (
541                 cloudI,
542                 lagrangianObjects,
543                 lagrangianVectorFields
544             );
546             lagrangianFieldDecomposer::readFields
547             (
548                 cloudI,
549                 lagrangianObjects,
550                 lagrangianSphericalTensorFields
551             );
553             lagrangianFieldDecomposer::readFields
554             (
555                 cloudI,
556                 lagrangianObjects,
557                 lagrangianSymmTensorFields
558             );
560             lagrangianFieldDecomposer::readFields
561             (
562                 cloudI,
563                 lagrangianObjects,
564                 lagrangianTensorFields
565             );
567             cloudI++;
568         }
569     }
571     lagrangianPositions.setSize(cloudI);
572     cellParticles.setSize(cloudI);
573     lagrangianLabelFields.setSize(cloudI);
574     lagrangianScalarFields.setSize(cloudI);
575     lagrangianVectorFields.setSize(cloudI);
576     lagrangianSphericalTensorFields.setSize(cloudI);
577     lagrangianSymmTensorFields.setSize(cloudI);
578     lagrangianTensorFields.setSize(cloudI);
581     // Any uniform data to copy/link?
582     fileName uniformDir("uniform");
584     if (isDir(runTime.timePath()/uniformDir))
585     {
586         Info<< "Detected additional non-decomposed files in "
587             << runTime.timePath()/uniformDir
588             << endl;
589     }
590     else
591     {
592         uniformDir.clear();
593     }
595     Info<< endl;
597     // Split the fields over processors
598     for (label procI = 0; procI < mesh.nProcs(); procI++)
599     {
600         Info<< "Processor " << procI << ": field transfer" << endl;
602         // open the database
603         Time processorDb
604         (
605             Time::controlDictName,
606             args.rootPath(),
607             args.caseName()/fileName(word("processor") + name(procI))
608         );
610         processorDb.setTime(runTime);
612         // Remove files remnants that can cause horrible problems
613         // - mut and nut are used to mark the new turbulence models,
614         //   their existence prevents old models from being upgraded
615         // 1.6.x merge.  HJ, 25/Aug/2010
616         {
617             fileName timeDir(processorDb.path()/processorDb.timeName());
619             rm(timeDir/"mut");
620             rm(timeDir/"nut");
621             rm(timeDir/"mut.gz");
622             rm(timeDir/"nut.gz");
623         }
625         // read the mesh
626         fvMesh procMesh
627         (
628             IOobject
629             (
630                 regionName,
631                 processorDb.timeName(),
632                 processorDb
633             )
634         );
636         labelIOList cellProcAddressing
637         (
638             IOobject
639             (
640                 "cellProcAddressing",
641                 procMesh.facesInstance(),
642                 procMesh.meshSubDir,
643                 procMesh,
644                 IOobject::MUST_READ,
645                 IOobject::NO_WRITE
646             )
647         );
649         labelIOList boundaryProcAddressing
650         (
651             IOobject
652             (
653                 "boundaryProcAddressing",
654                 procMesh.facesInstance(),
655                 procMesh.meshSubDir,
656                 procMesh,
657                 IOobject::MUST_READ,
658                 IOobject::NO_WRITE
659             )
660         );
662         // FV fields
663         if
664         (
665             volScalarFields.size()
666          || volVectorFields.size()
667          || volSphericalTensorFields.size()
668          || volSymmTensorFields.size()
669          || volTensorFields.size()
670          || surfaceScalarFields.size()
671          || surfaceVectorFields.size()
672          || surfaceSphericalTensorFields.size()
673          || surfaceSymmTensorFields.size()
674          || surfaceTensorFields.size()
675         )
676         {
677             labelIOList faceProcAddressing
678             (
679                 IOobject
680                 (
681                     "faceProcAddressing",
682                     procMesh.facesInstance(),
683                     procMesh.meshSubDir,
684                     procMesh,
685                     IOobject::MUST_READ,
686                     IOobject::NO_WRITE
687                 )
688             );
690             fvFieldDecomposer fieldDecomposer
691             (
692                 mesh,
693                 procMesh,
694                 faceProcAddressing,
695                 cellProcAddressing,
696                 boundaryProcAddressing
697             );
699             fieldDecomposer.decomposeFields(volScalarFields);
700             fieldDecomposer.decomposeFields(volVectorFields);
701             fieldDecomposer.decomposeFields(volSphericalTensorFields);
702             fieldDecomposer.decomposeFields(volSymmTensorFields);
703             fieldDecomposer.decomposeFields(volTensorFields);
705             fieldDecomposer.decomposeFields(surfaceScalarFields);
706             fieldDecomposer.decomposeFields(surfaceVectorFields);
707             fieldDecomposer.decomposeFields(surfaceSphericalTensorFields);
708             fieldDecomposer.decomposeFields(surfaceSymmTensorFields);
709             fieldDecomposer.decomposeFields(surfaceTensorFields);
710         }
713         // Point fields
714         if
715         (
716             pointScalarFields.size()
717          || pointVectorFields.size()
718          || pointSphericalTensorFields.size()
719          || pointSymmTensorFields.size()
720          || pointTensorFields.size()
721         )
722         {
723             labelIOList pointProcAddressing
724             (
725                 IOobject
726                 (
727                     "pointProcAddressing",
728                     procMesh.facesInstance(),
729                     procMesh.meshSubDir,
730                     procMesh,
731                     IOobject::MUST_READ,
732                     IOobject::NO_WRITE
733                 )
734             );
736             pointMesh procPMesh(procMesh, true);
738             pointFieldDecomposer fieldDecomposer
739             (
740                 pMesh,
741                 procPMesh,
742                 pointProcAddressing,
743                 boundaryProcAddressing
744             );
746             fieldDecomposer.decomposeFields(pointScalarFields);
747             fieldDecomposer.decomposeFields(pointVectorFields);
748             fieldDecomposer.decomposeFields(pointSphericalTensorFields);
749             fieldDecomposer.decomposeFields(pointSymmTensorFields);
750             fieldDecomposer.decomposeFields(pointTensorFields);
751         }
754         // tetPoint fields
755         if (tetMeshPtr)
756         {
757             const tetPolyMesh& tetMesh = *tetMeshPtr;
758             tetPolyMesh procTetMesh(procMesh);
760             // Read the point addressing information
761             labelIOList pointProcAddressing
762             (
763                 IOobject
764                 (
765                     "pointProcAddressing",
766                     procMesh.facesInstance(),
767                     procMesh.meshSubDir,
768                     procMesh,
769                     IOobject::MUST_READ,
770                     IOobject::NO_WRITE
771                 )
772             );
774             // Read the point addressing information
775             labelIOList faceProcAddressing
776             (
777                 IOobject
778                 (
779                     "faceProcAddressing",
780                     procMesh.facesInstance(),
781                     procMesh.meshSubDir,
782                     procMesh,
783                     IOobject::MUST_READ,
784                     IOobject::NO_WRITE
785                 )
786             );
788             tetPointFieldDecomposer fieldDecomposer
789             (
790                 tetMesh,
791                 procTetMesh,
792                 pointProcAddressing,
793                 faceProcAddressing,
794                 cellProcAddressing,
795                 boundaryProcAddressing
796             );
798             fieldDecomposer.decomposeFields(tetPointScalarFields);
799             fieldDecomposer.decomposeFields(tetPointVectorFields);
800             fieldDecomposer.decomposeFields(tetPointSphericalTensorFields);
801             fieldDecomposer.decomposeFields(tetPointSymmTensorFields);
802             fieldDecomposer.decomposeFields(tetPointTensorFields);
804             fieldDecomposer.decomposeFields(elementScalarFields);
805             fieldDecomposer.decomposeFields(elementVectorFields);
806         }
809         // If there is lagrangian data write it out
810         forAll(lagrangianPositions, cloudI)
811         {
812             if (lagrangianPositions[cloudI].size())
813             {
814                 lagrangianFieldDecomposer fieldDecomposer
815                 (
816                     mesh,
817                     procMesh,
818                     cellProcAddressing,
819                     cloudDirs[cloudI],
820                     lagrangianPositions[cloudI],
821                     cellParticles[cloudI]
822                 );
824                 // Lagrangian fields
825                 if
826                 (
827                     lagrangianLabelFields[cloudI].size()
828                  || lagrangianScalarFields[cloudI].size()
829                  || lagrangianVectorFields[cloudI].size()
830                  || lagrangianSphericalTensorFields[cloudI].size()
831                  || lagrangianSymmTensorFields[cloudI].size()
832                  || lagrangianTensorFields[cloudI].size()
833                 )
834                 {
835                     fieldDecomposer.decomposeFields
836                     (
837                         cloudDirs[cloudI],
838                         lagrangianLabelFields[cloudI]
839                     );
840                     fieldDecomposer.decomposeFields
841                     (
842                         cloudDirs[cloudI],
843                         lagrangianScalarFields[cloudI]
844                     );
845                     fieldDecomposer.decomposeFields
846                     (
847                         cloudDirs[cloudI],
848                         lagrangianVectorFields[cloudI]
849                     );
850                     fieldDecomposer.decomposeFields
851                     (
852                         cloudDirs[cloudI],
853                         lagrangianSphericalTensorFields[cloudI]
854                     );
855                     fieldDecomposer.decomposeFields
856                     (
857                         cloudDirs[cloudI],
858                         lagrangianSymmTensorFields[cloudI]
859                     );
860                     fieldDecomposer.decomposeFields
861                     (
862                         cloudDirs[cloudI],
863                         lagrangianTensorFields[cloudI]
864                     );
865                 }
866             }
867         }
870         // Any non-decomposed data to copy?
871         if (uniformDir.size())
872         {
873             const fileName timePath = processorDb.timePath();
875             if (copyUniform || mesh.distributed())
876             {
877                 cp
878                 (
879                     runTime.timePath()/uniformDir,
880                     timePath/uniformDir
881                 );
882             }
883             else
884             {
885                 // Link with relative paths
886                 const string parentPath = string("..")/"..";
888                 fileName currentDir(cwd());
889                 chDir(timePath);
890                 if (!exists(uniformDir))
891                 {
892                     ln
893                     (
894                         parentPath/runTime.timeName()/uniformDir,
895                         uniformDir
896                     );
897                     chDir(currentDir);
898                 }
899             }
900         }
901     }
904     if (tetMeshPtr)
905     {
906         delete tetMeshPtr;
907         tetMeshPtr = NULL;
908     }
911     // Finite area mesh and field decomposition
913     IOobject faMeshBoundaryIOobj
914     (
915         "faBoundary",
916         mesh.time().findInstance
917         (
918             mesh.dbDir()/polyMesh::meshSubDir,
919             "boundary"
920         ),
921         faMesh::meshSubDir,
922         mesh,
923         IOobject::READ_IF_PRESENT,
924         IOobject::NO_WRITE
925     );
928     if(faMeshBoundaryIOobj.headerOk())
929     {
930         Info << "\nFinite area mesh decomposition" << endl;
932         faMeshDecomposition aMesh(mesh);
934         aMesh.decomposeMesh(filterPatches);
936         aMesh.writeDecomposition();
939         // Construct the area fields
940         // ~~~~~~~~~~~~~~~~~~~~~~~~
941         PtrList<areaScalarField> areaScalarFields;
942         readFields(aMesh, objects, areaScalarFields);
944         PtrList<areaVectorField> areaVectorFields;
945         readFields(aMesh, objects, areaVectorFields);
947         PtrList<areaSphericalTensorField> areaSphericalTensorFields;
948         readFields(aMesh, objects, areaSphericalTensorFields);
950         PtrList<areaSymmTensorField> areaSymmTensorFields;
951         readFields(aMesh, objects, areaSymmTensorFields);
953         PtrList<areaTensorField> areaTensorFields;
954         readFields(aMesh, objects, areaTensorFields);
957         // Construct the edge fields
958         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
959         PtrList<edgeScalarField> edgeScalarFields;
960         readFields(aMesh, objects, edgeScalarFields);
962         Info << endl;
964         // Split the fields over processors
965         for (label procI = 0; procI < mesh.nProcs(); procI++)
966         {
967             Info<< "Processor " << procI
968                 << ": finite area field transfer" << endl;
970             // open the database
971             Time processorDb
972             (
973                 Time::controlDictName,
974                 args.rootPath(),
975                 args.caseName()/fileName(word("processor") + name(procI))
976             );
978             processorDb.setTime(runTime);
980             // Read the mesh
981             fvMesh procFvMesh
982             (
983                 IOobject
984                 (
985                     regionName,
986                     processorDb.timeName(),
987                     processorDb
988                 )
989             );
991             faMesh procMesh(procFvMesh);
993             labelIOList faceProcAddressing
994             (
995                 IOobject
996                 (
997                     "faceProcAddressing",
998                     "constant",
999                     procMesh.meshSubDir,
1000                     procFvMesh,
1001                     IOobject::MUST_READ,
1002                     IOobject::NO_WRITE
1003                 )
1004             );
1006             labelIOList boundaryProcAddressing
1007             (
1008                 IOobject
1009                 (
1010                     "boundaryProcAddressing",
1011                     "constant",
1012                     procMesh.meshSubDir,
1013                     procFvMesh,
1014                     IOobject::MUST_READ,
1015                     IOobject::NO_WRITE
1016                 )
1017             );
1019             // FA fields
1020             if
1021             (
1022                 areaScalarFields.size()
1023              || areaVectorFields.size()
1024              || areaSphericalTensorFields.size()
1025              || areaSymmTensorFields.size()
1026              || areaTensorFields.size()
1027              || edgeScalarFields.size()
1028             )
1029             {
1030                 labelIOList edgeProcAddressing
1031                 (
1032                     IOobject
1033                     (
1034                         "edgeProcAddressing",
1035                         "constant",
1036                         procMesh.meshSubDir,
1037                         procFvMesh,
1038                         IOobject::MUST_READ,
1039                         IOobject::NO_WRITE
1040                     )
1041                 );
1043                 faFieldDecomposer fieldDecomposer
1044                 (
1045                     aMesh,
1046                     procMesh,
1047                     edgeProcAddressing,
1048                     faceProcAddressing,
1049                     boundaryProcAddressing
1050                 );
1052                 fieldDecomposer.decomposeFields(areaScalarFields);
1053                 fieldDecomposer.decomposeFields(areaVectorFields);
1054                 fieldDecomposer.decomposeFields(areaSphericalTensorFields);
1055                 fieldDecomposer.decomposeFields(areaSymmTensorFields);
1056                 fieldDecomposer.decomposeFields(areaTensorFields);
1058                 fieldDecomposer.decomposeFields(edgeScalarFields);
1059             }
1060         }
1061     }
1064     Info<< "\nEnd.\n" << endl;
1066     return 0;
1070 // ************************************************************************* //