ENH: decomposePar.C: add shortcircuit to avoid allocating point mappers
[OpenFOAM-2.0.x.git] / applications / utilities / parallelProcessing / decomposePar / decomposePar.C
blob02acd59236545782ba806d64b01ae7c8ae2496e2
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     decomposePar
27 Description
28     Automatically decomposes a mesh and fields of a case for parallel
29     execution of OpenFOAM.
31 Usage
33     - decomposePar [OPTION]
35     \param -cellDist \n
36     Write the cell distribution as a labelList, for use with 'manual'
37     decomposition method or as a volScalarField for post-processing.
39     \param -region regionName \n
40     Decompose named region. Does not check for existence of processor*.
42     \param -copyUniform \n
43     Copy any \a uniform directories too.
45     \param -constant \n
46     Override controlDict settings and use constant directory.
48     \param -fields \n
49     Use existing geometry decomposition and convert fields only.
51     \param -force \n
52     Remove any existing \a processor subdirectories before decomposing the
53     geometry.
55     \param -ifRequired \n
56     Only decompose the geometry if the number of domains has changed from a
57     previous decomposition. No \a processor subdirectories will be removed
58     unless the \a -force option is also specified. This option can be used
59     to avoid redundant geometry decomposition (eg, in scripts), but should
60     be used with caution when the underlying (serial) geometry or the
61     decomposition method etc. have been changed between decompositions.
63 \*---------------------------------------------------------------------------*/
65 #include "OSspecific.H"
66 #include "fvCFD.H"
67 #include "IOobjectList.H"
68 #include "domainDecomposition.H"
69 #include "labelIOField.H"
70 #include "labelFieldIOField.H"
71 #include "scalarIOField.H"
72 #include "scalarFieldIOField.H"
73 #include "vectorIOField.H"
74 #include "vectorFieldIOField.H"
75 #include "sphericalTensorIOField.H"
76 #include "sphericalTensorFieldIOField.H"
77 #include "symmTensorIOField.H"
78 #include "symmTensorFieldIOField.H"
79 #include "tensorIOField.H"
80 #include "tensorFieldIOField.H"
81 #include "pointFields.H"
83 #include "readFields.H"
84 #include "dimFieldDecomposer.H"
85 #include "fvFieldDecomposer.H"
86 #include "pointFieldDecomposer.H"
87 #include "lagrangianFieldDecomposer.H"
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 int main(int argc, char *argv[])
93     argList::addNote
94     (
95         "decompose a mesh and fields of a case for parallel execution"
96     );
98     argList::noParallel();
99     #include "addRegionOption.H"
100     argList::addBoolOption
101     (
102         "cellDist",
103         "write cell distribution as a labelList - for use with 'manual' "
104         "decomposition method or as a volScalarField for post-processing."
105     );
106     argList::addBoolOption
107     (
108         "copyUniform",
109         "copy any uniform/ directories too"
110     );
111     argList::addBoolOption
112     (
113         "fields",
114         "use existing geometry decomposition and convert fields only"
115     );
116     argList::addBoolOption
117     (
118         "force",
119         "remove existing processor*/ subdirs before decomposing the geometry"
120     );
121     argList::addBoolOption
122     (
123         "ifRequired",
124         "only decompose geometry if the number of domains has changed"
125     );
126     argList::addBoolOption
127     (
128         "constant",
129         "include the 'constant/' dir in the times list"
130     );
132     #include "setRootCase.H"
134     word regionName = fvMesh::defaultRegion;
135     word regionDir = word::null;
137     if (args.optionReadIfPresent("region", regionName))
138     {
139         regionDir = regionName;
140         Info<< "Decomposing mesh " << regionName << nl << endl;
141     }
143     bool writeCellDist           = args.optionFound("cellDist");
144     bool copyUniform             = args.optionFound("copyUniform");
145     bool decomposeFieldsOnly     = args.optionFound("fields");
146     bool forceOverwrite          = args.optionFound("force");
147     bool ifRequiredDecomposition = args.optionFound("ifRequired");
149     #include "createTime.H"
151     // Allow -constant to override controlDict settings.
152     if (args.optionFound("constant"))
153     {
154         instantList timeDirs = timeSelector::select0(runTime, args);
155         if (runTime.timeName() != runTime.constant())
156         {
157             FatalErrorIn(args.executable())
158                 << "No '" << runTime.constant() << "' time present." << endl
159                 << "Valid times are " << runTime.times()
160                 << exit(FatalError);
161         }
162     }
165     Info<< "Time = " << runTime.timeName() << endl;
167     // determine the existing processor count directly
168     label nProcs = 0;
169     while
170     (
171         isDir
172         (
173             runTime.path()
174           / (word("processor") + name(nProcs))
175           / runTime.constant()
176           / regionDir
177           / polyMesh::meshSubDir
178         )
179     )
180     {
181         ++nProcs;
182     }
184     // get requested numberOfSubdomains
185     const label nDomains = readLabel
186     (
187         IOdictionary
188         (
189             IOobject
190             (
191                 "decomposeParDict",
192                 runTime.time().system(),
193                 regionDir,          // use region if non-standard
194                 runTime,
195                 IOobject::MUST_READ_IF_MODIFIED,
196                 IOobject::NO_WRITE,
197                 false
198             )
199         ).lookup("numberOfSubdomains")
200     );
202     if (decomposeFieldsOnly)
203     {
204         // Sanity check on previously decomposed case
205         if (nProcs != nDomains)
206         {
207             FatalErrorIn(args.executable())
208                 << "Specified -fields, but the case was decomposed with "
209                 << nProcs << " domains"
210                 << nl
211                 << "instead of " << nDomains
212                 << " domains as specified in decomposeParDict"
213                 << nl
214                 << exit(FatalError);
215         }
216     }
217     else if (nProcs)
218     {
219         bool procDirsProblem = true;
221         if (ifRequiredDecomposition && nProcs == nDomains)
222         {
223             // we can reuse the decomposition
224             decomposeFieldsOnly = true;
225             procDirsProblem = false;
226             forceOverwrite = false;
228             Info<< "Using existing processor directories" << nl;
229         }
231         if (forceOverwrite)
232         {
233             Info<< "Removing " << nProcs
234                 << " existing processor directories" << endl;
236             // remove existing processor dirs
237             // reverse order to avoid gaps if someone interrupts the process
238             for (label procI = nProcs-1; procI >= 0; --procI)
239             {
240                 fileName procDir
241                 (
242                     runTime.path()/(word("processor") + name(procI))
243                 );
245                 rmDir(procDir);
246             }
248             procDirsProblem = false;
249         }
251         if (procDirsProblem)
252         {
253             FatalErrorIn(args.executable())
254                 << "Case is already decomposed with " << nProcs
255                 << " domains, use the -force option or manually" << nl
256                 << "remove processor directories before decomposing. e.g.,"
257                 << nl
258                 << "    rm -rf " << runTime.path().c_str() << "/processor*"
259                 << nl
260                 << exit(FatalError);
261         }
262     }
264     Info<< "Create mesh" << endl;
265     domainDecomposition mesh
266     (
267         IOobject
268         (
269             regionName,
270             runTime.timeName(),
271             runTime
272         )
273     );
275     // Decompose the mesh
276     if (!decomposeFieldsOnly)
277     {
278         mesh.decomposeMesh();
280         mesh.writeDecomposition();
282         if (writeCellDist)
283         {
284             const labelList& procIds = mesh.cellToProc();
286             // Write the decomposition as labelList for use with 'manual'
287             // decomposition method.
288             labelIOList cellDecomposition
289             (
290                 IOobject
291                 (
292                     "cellDecomposition",
293                     mesh.facesInstance(),
294                     mesh,
295                     IOobject::NO_READ,
296                     IOobject::NO_WRITE,
297                     false
298                 ),
299                 procIds
300             );
301             cellDecomposition.write();
303             Info<< nl << "Wrote decomposition to "
304                 << cellDecomposition.objectPath()
305                 << " for use in manual decomposition." << endl;
307             // Write as volScalarField for postprocessing.
308             volScalarField cellDist
309             (
310                 IOobject
311                 (
312                     "cellDist",
313                     runTime.timeName(),
314                     mesh,
315                     IOobject::NO_READ,
316                     IOobject::AUTO_WRITE
317                 ),
318                 mesh,
319                 dimensionedScalar("cellDist", dimless, 0),
320                 zeroGradientFvPatchScalarField::typeName
321             );
323             forAll(procIds, celli)
324             {
325                cellDist[celli] = procIds[celli];
326             }
328             cellDist.write();
330             Info<< nl << "Wrote decomposition as volScalarField to "
331                 << cellDist.name() << " for use in postprocessing."
332                 << endl;
333         }
334     }
337     // Search for list of objects for this time
338     IOobjectList objects(mesh, runTime.timeName());
340     // Construct the vol fields
341     // ~~~~~~~~~~~~~~~~~~~~~~~~
342     PtrList<volScalarField> volScalarFields;
343     readFields(mesh, objects, volScalarFields);
345     PtrList<volVectorField> volVectorFields;
346     readFields(mesh, objects, volVectorFields);
348     PtrList<volSphericalTensorField> volSphericalTensorFields;
349     readFields(mesh, objects, volSphericalTensorFields);
351     PtrList<volSymmTensorField> volSymmTensorFields;
352     readFields(mesh, objects, volSymmTensorFields);
354     PtrList<volTensorField> volTensorFields;
355     readFields(mesh, objects, volTensorFields);
358     // Construct the dimensioned fields
359     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
360     PtrList<DimensionedField<scalar, volMesh> > dimScalarFields;
361     readFields(mesh, objects, dimScalarFields);
363     PtrList<DimensionedField<vector, volMesh> > dimVectorFields;
364     readFields(mesh, objects, dimVectorFields);
366     PtrList<DimensionedField<sphericalTensor, volMesh> >
367         dimSphericalTensorFields;
368     readFields(mesh, objects, dimSphericalTensorFields);
370     PtrList<DimensionedField<symmTensor, volMesh> > dimSymmTensorFields;
371     readFields(mesh, objects, dimSymmTensorFields);
373     PtrList<DimensionedField<tensor, volMesh> > dimTensorFields;
374     readFields(mesh, objects, dimTensorFields);
377     // Construct the surface fields
378     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
379     PtrList<surfaceScalarField> surfaceScalarFields;
380     readFields(mesh, objects, surfaceScalarFields);
381     PtrList<surfaceVectorField> surfaceVectorFields;
382     readFields(mesh, objects, surfaceVectorFields);
383     PtrList<surfaceSphericalTensorField> surfaceSphericalTensorFields;
384     readFields(mesh, objects, surfaceSphericalTensorFields);
385     PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
386     readFields(mesh, objects, surfaceSymmTensorFields);
387     PtrList<surfaceTensorField> surfaceTensorFields;
388     readFields(mesh, objects, surfaceTensorFields);
391     // Construct the point fields
392     // ~~~~~~~~~~~~~~~~~~~~~~~~~~
393     pointMesh pMesh(mesh);
395     PtrList<pointScalarField> pointScalarFields;
396     readFields(pMesh, objects, pointScalarFields);
398     PtrList<pointVectorField> pointVectorFields;
399     readFields(pMesh, objects, pointVectorFields);
401     PtrList<pointSphericalTensorField> pointSphericalTensorFields;
402     readFields(pMesh, objects, pointSphericalTensorFields);
404     PtrList<pointSymmTensorField> pointSymmTensorFields;
405     readFields(pMesh, objects, pointSymmTensorFields);
407     PtrList<pointTensorField> pointTensorFields;
408     readFields(pMesh, objects, pointTensorFields);
411     // Construct the Lagrangian fields
412     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
414     fileNameList cloudDirs
415     (
416         readDir(runTime.timePath()/cloud::prefix, fileName::DIRECTORY)
417     );
419     // Particles
420     PtrList<Cloud<indexedParticle> > lagrangianPositions(cloudDirs.size());
421     // Particles per cell
422     PtrList< List<SLList<indexedParticle*>*> > cellParticles(cloudDirs.size());
424     PtrList<PtrList<labelIOField> > lagrangianLabelFields(cloudDirs.size());
425     PtrList<PtrList<labelFieldCompactIOField> > lagrangianLabelFieldFields
426     (
427         cloudDirs.size()
428     );
429     PtrList<PtrList<scalarIOField> > lagrangianScalarFields(cloudDirs.size());
430     PtrList<PtrList<scalarFieldCompactIOField> > lagrangianScalarFieldFields
431     (
432         cloudDirs.size()
433     );
434     PtrList<PtrList<vectorIOField> > lagrangianVectorFields(cloudDirs.size());
435     PtrList<PtrList<vectorFieldCompactIOField> > lagrangianVectorFieldFields
436     (
437         cloudDirs.size()
438     );
439     PtrList<PtrList<sphericalTensorIOField> > lagrangianSphericalTensorFields
440     (
441         cloudDirs.size()
442     );
443     PtrList<PtrList<sphericalTensorFieldCompactIOField> >
444         lagrangianSphericalTensorFieldFields(cloudDirs.size());
445     PtrList<PtrList<symmTensorIOField> > lagrangianSymmTensorFields
446     (
447         cloudDirs.size()
448     );
449     PtrList<PtrList<symmTensorFieldCompactIOField> >
450     lagrangianSymmTensorFieldFields
451     (
452         cloudDirs.size()
453     );
454     PtrList<PtrList<tensorIOField> > lagrangianTensorFields
455     (
456         cloudDirs.size()
457     );
458     PtrList<PtrList<tensorFieldCompactIOField> > lagrangianTensorFieldFields
459     (
460         cloudDirs.size()
461     );
463     label cloudI = 0;
465     forAll(cloudDirs, i)
466     {
467         IOobjectList sprayObjs
468         (
469             mesh,
470             runTime.timeName(),
471             cloud::prefix/cloudDirs[i]
472         );
474         IOobject* positionsPtr = sprayObjs.lookup("positions");
476         if (positionsPtr)
477         {
478             // Read lagrangian particles
479             // ~~~~~~~~~~~~~~~~~~~~~~~~~
481             Info<< "Identified lagrangian data set: " << cloudDirs[i] << endl;
483             lagrangianPositions.set
484             (
485                 cloudI,
486                 new Cloud<indexedParticle>
487                 (
488                     mesh,
489                     cloudDirs[i],
490                     false
491                 )
492             );
495             // Sort particles per cell
496             // ~~~~~~~~~~~~~~~~~~~~~~~
498             cellParticles.set
499             (
500                 cloudI,
501                 new List<SLList<indexedParticle*>*>
502                 (
503                     mesh.nCells(),
504                     static_cast<SLList<indexedParticle*>*>(NULL)
505                 )
506             );
508             label i = 0;
510             forAllIter
511             (
512                 Cloud<indexedParticle>,
513                 lagrangianPositions[cloudI],
514                 iter
515             )
516             {
517                 iter().index() = i++;
519                 label celli = iter().cell();
521                 // Check
522                 if (celli < 0 || celli >= mesh.nCells())
523                 {
524                     FatalErrorIn(args.executable())
525                         << "Illegal cell number " << celli
526                         << " for particle with index " << iter().index()
527                         << " at position " << iter().position() << nl
528                         << "Cell number should be between 0 and "
529                         << mesh.nCells()-1 << nl
530                         << "On this mesh the particle should be in cell "
531                         << mesh.findCell(iter().position())
532                         << exit(FatalError);
533                 }
535                 if (!cellParticles[cloudI][celli])
536                 {
537                     cellParticles[cloudI][celli] = new SLList<indexedParticle*>
538                     ();
539                 }
541                 cellParticles[cloudI][celli]->append(&iter());
542             }
544             // Read fields
545             // ~~~~~~~~~~~
547             IOobjectList lagrangianObjects
548             (
549                 mesh,
550                 runTime.timeName(),
551                 cloud::prefix/cloudDirs[cloudI]
552             );
554             lagrangianFieldDecomposer::readFields
555             (
556                 cloudI,
557                 lagrangianObjects,
558                 lagrangianLabelFields
559             );
561             lagrangianFieldDecomposer::readFieldFields
562             (
563                 cloudI,
564                 lagrangianObjects,
565                 lagrangianLabelFieldFields
566             );
568             lagrangianFieldDecomposer::readFields
569             (
570                 cloudI,
571                 lagrangianObjects,
572                 lagrangianScalarFields
573             );
575             lagrangianFieldDecomposer::readFieldFields
576             (
577                 cloudI,
578                 lagrangianObjects,
579                 lagrangianScalarFieldFields
580             );
582             lagrangianFieldDecomposer::readFields
583             (
584                 cloudI,
585                 lagrangianObjects,
586                 lagrangianVectorFields
587             );
589             lagrangianFieldDecomposer::readFieldFields
590             (
591                 cloudI,
592                 lagrangianObjects,
593                 lagrangianVectorFieldFields
594             );
596             lagrangianFieldDecomposer::readFields
597             (
598                 cloudI,
599                 lagrangianObjects,
600                 lagrangianSphericalTensorFields
601             );
603             lagrangianFieldDecomposer::readFieldFields
604             (
605                 cloudI,
606                 lagrangianObjects,
607                 lagrangianSphericalTensorFieldFields
608             );
610             lagrangianFieldDecomposer::readFields
611             (
612                 cloudI,
613                 lagrangianObjects,
614                 lagrangianSymmTensorFields
615             );
617             lagrangianFieldDecomposer::readFieldFields
618             (
619                 cloudI,
620                 lagrangianObjects,
621                 lagrangianSymmTensorFieldFields
622             );
624             lagrangianFieldDecomposer::readFields
625             (
626                 cloudI,
627                 lagrangianObjects,
628                 lagrangianTensorFields
629             );
631             lagrangianFieldDecomposer::readFieldFields
632             (
633                 cloudI,
634                 lagrangianObjects,
635                 lagrangianTensorFieldFields
636             );
638             cloudI++;
639         }
640     }
642     lagrangianPositions.setSize(cloudI);
643     cellParticles.setSize(cloudI);
644     lagrangianLabelFields.setSize(cloudI);
645     lagrangianLabelFieldFields.setSize(cloudI);
646     lagrangianScalarFields.setSize(cloudI);
647     lagrangianScalarFieldFields.setSize(cloudI);
648     lagrangianVectorFields.setSize(cloudI);
649     lagrangianVectorFieldFields.setSize(cloudI);
650     lagrangianSphericalTensorFields.setSize(cloudI);
651     lagrangianSphericalTensorFieldFields.setSize(cloudI);
652     lagrangianSymmTensorFields.setSize(cloudI);
653     lagrangianSymmTensorFieldFields.setSize(cloudI);
654     lagrangianTensorFields.setSize(cloudI);
655     lagrangianTensorFieldFields.setSize(cloudI);
658     // Any uniform data to copy/link?
659     fileName uniformDir("uniform");
661     if (isDir(runTime.timePath()/uniformDir))
662     {
663         Info<< "Detected additional non-decomposed files in "
664             << runTime.timePath()/uniformDir
665             << endl;
666     }
667     else
668     {
669         uniformDir.clear();
670     }
672     Info<< endl;
674     // split the fields over processors
675     for (label procI = 0; procI < mesh.nProcs(); procI++)
676     {
677         Info<< "Processor " << procI << ": field transfer" << endl;
679         // open the database
680         Time processorDb
681         (
682             Time::controlDictName,
683             args.rootPath(),
684             args.caseName()/fileName(word("processor") + name(procI))
685         );
687         processorDb.setTime(runTime);
689         // remove files remnants that can cause horrible problems
690         // - mut and nut are used to mark the new turbulence models,
691         //   their existence prevents old models from being upgraded
692         {
693             fileName timeDir(processorDb.path()/processorDb.timeName());
695             rm(timeDir/"mut");
696             rm(timeDir/"nut");
697         }
699         // read the mesh
700         fvMesh procMesh
701         (
702             IOobject
703             (
704                 regionName,
705                 processorDb.timeName(),
706                 processorDb
707             )
708         );
710         labelIOList faceProcAddressing
711         (
712             IOobject
713             (
714                 "faceProcAddressing",
715                 procMesh.facesInstance(),
716                 procMesh.meshSubDir,
717                 procMesh,
718                 IOobject::MUST_READ,
719                 IOobject::NO_WRITE
720             )
721         );
723         labelIOList cellProcAddressing
724         (
725             IOobject
726             (
727                 "cellProcAddressing",
728                 procMesh.facesInstance(),
729                 procMesh.meshSubDir,
730                 procMesh,
731                 IOobject::MUST_READ,
732                 IOobject::NO_WRITE
733             )
734         );
736         labelIOList boundaryProcAddressing
737         (
738             IOobject
739             (
740                 "boundaryProcAddressing",
741                 procMesh.facesInstance(),
742                 procMesh.meshSubDir,
743                 procMesh,
744                 IOobject::MUST_READ,
745                 IOobject::NO_WRITE
746             )
747         );
749         // FV fields
750         {
751             fvFieldDecomposer fieldDecomposer
752             (
753                 mesh,
754                 procMesh,
755                 faceProcAddressing,
756                 cellProcAddressing,
757                 boundaryProcAddressing
758             );
760             fieldDecomposer.decomposeFields(volScalarFields);
761             fieldDecomposer.decomposeFields(volVectorFields);
762             fieldDecomposer.decomposeFields(volSphericalTensorFields);
763             fieldDecomposer.decomposeFields(volSymmTensorFields);
764             fieldDecomposer.decomposeFields(volTensorFields);
766             fieldDecomposer.decomposeFields(surfaceScalarFields);
767             fieldDecomposer.decomposeFields(surfaceVectorFields);
768             fieldDecomposer.decomposeFields(surfaceSphericalTensorFields);
769             fieldDecomposer.decomposeFields(surfaceSymmTensorFields);
770             fieldDecomposer.decomposeFields(surfaceTensorFields);
771         }
773         // Dimensioned fields
774         {
775             dimFieldDecomposer fieldDecomposer
776             (
777                 mesh,
778                 procMesh,
779                 faceProcAddressing,
780                 cellProcAddressing
781             );
783             fieldDecomposer.decomposeFields(dimScalarFields);
784             fieldDecomposer.decomposeFields(dimVectorFields);
785             fieldDecomposer.decomposeFields(dimSphericalTensorFields);
786             fieldDecomposer.decomposeFields(dimSymmTensorFields);
787             fieldDecomposer.decomposeFields(dimTensorFields);
788         }
791         // Point fields
792         if
793         (
794             pointScalarFields.size()
795          || pointVectorFields.size()
796          || pointSphericalTensorFields.size()
797          || pointSymmTensorFields.size()
798          || pointTensorFields.size()
799         )
800         {
801             labelIOList pointProcAddressing
802             (
803                 IOobject
804                 (
805                     "pointProcAddressing",
806                     procMesh.facesInstance(),
807                     procMesh.meshSubDir,
808                     procMesh,
809                     IOobject::MUST_READ,
810                     IOobject::NO_WRITE
811                 )
812             );
814             pointMesh procPMesh(procMesh);
816             pointFieldDecomposer fieldDecomposer
817             (
818                 pMesh,
819                 procPMesh,
820                 pointProcAddressing,
821                 boundaryProcAddressing
822             );
824             fieldDecomposer.decomposeFields(pointScalarFields);
825             fieldDecomposer.decomposeFields(pointVectorFields);
826             fieldDecomposer.decomposeFields(pointSphericalTensorFields);
827             fieldDecomposer.decomposeFields(pointSymmTensorFields);
828             fieldDecomposer.decomposeFields(pointTensorFields);
829         }
832         // If there is lagrangian data write it out
833         forAll(lagrangianPositions, cloudI)
834         {
835             if (lagrangianPositions[cloudI].size())
836             {
837                 lagrangianFieldDecomposer fieldDecomposer
838                 (
839                     mesh,
840                     procMesh,
841                     faceProcAddressing,
842                     cellProcAddressing,
843                     cloudDirs[cloudI],
844                     lagrangianPositions[cloudI],
845                     cellParticles[cloudI]
846                 );
848                 // Lagrangian fields
849                 {
850                     fieldDecomposer.decomposeFields
851                     (
852                         cloudDirs[cloudI],
853                         lagrangianLabelFields[cloudI]
854                     );
855                     fieldDecomposer.decomposeFieldFields
856                     (
857                         cloudDirs[cloudI],
858                         lagrangianLabelFieldFields[cloudI]
859                     );
860                     fieldDecomposer.decomposeFields
861                     (
862                         cloudDirs[cloudI],
863                         lagrangianScalarFields[cloudI]
864                     );
865                     fieldDecomposer.decomposeFieldFields
866                     (
867                         cloudDirs[cloudI],
868                         lagrangianScalarFieldFields[cloudI]
869                     );
870                     fieldDecomposer.decomposeFields
871                     (
872                         cloudDirs[cloudI],
873                         lagrangianVectorFields[cloudI]
874                     );
875                     fieldDecomposer.decomposeFieldFields
876                     (
877                         cloudDirs[cloudI],
878                         lagrangianVectorFieldFields[cloudI]
879                     );
880                     fieldDecomposer.decomposeFields
881                     (
882                         cloudDirs[cloudI],
883                         lagrangianSphericalTensorFields[cloudI]
884                     );
885                     fieldDecomposer.decomposeFieldFields
886                     (
887                         cloudDirs[cloudI],
888                         lagrangianSphericalTensorFieldFields[cloudI]
889                     );
890                     fieldDecomposer.decomposeFields
891                     (
892                         cloudDirs[cloudI],
893                         lagrangianSymmTensorFields[cloudI]
894                     );
895                     fieldDecomposer.decomposeFieldFields
896                     (
897                         cloudDirs[cloudI],
898                         lagrangianSymmTensorFieldFields[cloudI]
899                     );
900                     fieldDecomposer.decomposeFields
901                     (
902                         cloudDirs[cloudI],
903                         lagrangianTensorFields[cloudI]
904                     );
905                     fieldDecomposer.decomposeFieldFields
906                     (
907                         cloudDirs[cloudI],
908                         lagrangianTensorFieldFields[cloudI]
909                     );
910                 }
911             }
912         }
915         // Any non-decomposed data to copy?
916         if (uniformDir.size())
917         {
918             const fileName timePath = processorDb.timePath();
920             if (copyUniform || mesh.distributed())
921             {
922                 cp
923                 (
924                     runTime.timePath()/uniformDir,
925                     timePath/uniformDir
926                 );
927             }
928             else
929             {
930                 // link with relative paths
931                 const string parentPath = string("..")/"..";
933                 fileName currentDir(cwd());
934                 chDir(timePath);
935                 ln
936                 (
937                     parentPath/runTime.timeName()/uniformDir,
938                     uniformDir
939                 );
940                 chDir(currentDir);
941             }
942         }
943     }
946     Info<< "\nEnd.\n" << endl;
948     return 0;
952 // ************************************************************************* //