Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / applications / utilities / parallelProcessing / decomposePar / decomposePar.C
blob6fe3466e1cde05e440158a943fc30d3040ea4b3b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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 FOAM.
31 Usage
33     - decomposePar [OPTION]
35     @param -cellDist \n
36     Write the cell distribution as a labelList for use with 'manual'
37     decomposition method and 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 -fields \n
46     Use existing geometry decomposition and convert fields only.
48     @param -filterPatches \n
49     Remove empty patches when decomposing the geometry.
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 "processorFvPatchFields.H"
69 #include "domainDecomposition.H"
70 #include "labelIOField.H"
71 #include "scalarIOField.H"
72 #include "vectorIOField.H"
73 #include "sphericalTensorIOField.H"
74 #include "symmTensorIOField.H"
75 #include "tensorIOField.H"
77 #include "tetPointFields.H"
78 #include "elementFields.H"
79 #include "tetFemMatrices.H"
80 #include "tetPointFieldDecomposer.H"
82 #include "pointFields.H"
84 #include "readFields.H"
85 #include "fvFieldDecomposer.H"
86 #include "pointFieldDecomposer.H"
87 #include "lagrangianFieldDecomposer.H"
89 #include "faCFD.H"
90 #include "emptyFaPatch.H"
91 #include "faMeshDecomposition.H"
92 #include "faFieldDecomposer.H"
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
96 int main(int argc, char *argv[])
98     argList::noParallel();
99 #   include "addRegionOption.H"
100     argList::validOptions.insert("cellDist", "");
101     argList::validOptions.insert("copyUniform", "");
102     argList::validOptions.insert("fields", "");
103     argList::validOptions.insert("filterPatches", "");
104     argList::validOptions.insert("force", "");
105     argList::validOptions.insert("ifRequired", "");
107 #   include "setRootCase.H"
109     word regionName = fvMesh::defaultRegion;
110     word regionDir = word::null;
112     if (args.optionFound("region"))
113     {
114         regionName = args.option("region");
115         regionDir = regionName;
116         Info<< "Decomposing mesh " << regionName << nl << endl;
117     }
120     bool writeCellDist = args.optionFound("cellDist");
121     bool copyUniform = args.optionFound("copyUniform");
122     bool decomposeFieldsOnly = args.optionFound("fields");
123     bool filterPatches = args.optionFound("filterPatches");
124     bool forceOverwrite = args.optionFound("force");
125     bool ifRequiredDecomposition = args.optionFound("ifRequired");
127 #   include "createTime.H"
129     Info<< "Time = " << runTime.timeName() << endl;
131     // Determine the existing processor count directly
132     label nProcs = 0;
133     while
134     (
135         isDir
136         (
137             runTime.path()
138            /(word("processor") + name(nProcs))
139            /runTime.constant()
140            /regionDir
141            /polyMesh::meshSubDir
142         )
143     )
144     {
145         ++nProcs;
146     }
148     // get requested numberOfSubdomains
149     label nDomains = 0;
150     {
151         IOdictionary decompDict
152         (
153             IOobject
154             (
155                 "decomposeParDict",
156                 runTime.time().system(),
157                 regionDir,          // use region if non-standard
158                 runTime,
159                 IOobject::MUST_READ,
160                 IOobject::NO_WRITE,
161                 false
162             )
163         );
165         decompDict.lookup("numberOfSubdomains") >> nDomains;
166     }
168     if (decomposeFieldsOnly)
169     {
170         // Sanity check on previously decomposed case
171         if (nProcs != nDomains)
172         {
173             FatalErrorIn(args.executable())
174                 << "Specified -fields, but the case was decomposed with "
175                 << nProcs << " domains"
176                 << nl
177                 << "instead of " << nDomains
178                 << " domains as specified in decomposeParDict"
179                 << nl
180                 << exit(FatalError);
181         }
182     }
183     else if (nProcs)
184     {
185         bool procDirsProblem = true;
187         if (regionName != fvMesh::defaultRegion)
188         {
189             decomposeFieldsOnly = false;
190             procDirsProblem = false;
191         }
194         if (ifRequiredDecomposition && nProcs == nDomains)
195         {
196             // we can reuse the decomposition
197             decomposeFieldsOnly = true;
198             procDirsProblem = false;
199             forceOverwrite = false;
201             Info<< "Using existing processor directories" << nl;
202         }
204         if (forceOverwrite)
205         {
206             Info<< "Removing " << nProcs
207                 << " existing processor directories" << endl;
209             // remove existing processor dirs
210             // reverse order to avoid gaps if someone interrupts the process
211             for (label procI = nProcs-1; procI >= 0; --procI)
212             {
213                 fileName procDir
214                 (
215                     runTime.path()/(word("processor") + name(procI))
216                 );
218                 rmDir(procDir);
219             }
221             procDirsProblem = false;
222         }
224         if (procDirsProblem)
225         {
226             FatalErrorIn(args.executable())
227                 << "Case is already decomposed with " << nProcs
228                 << " domains, use the -force option or manually" << nl
229                 << "remove processor directories before decomposing. e.g.,"
230                 << nl
231                 << "    rm -rf " << runTime.path().c_str() << "/processor*"
232                 << nl
233                 << exit(FatalError);
234         }
235     }
237     Info<< "Create mesh for region " << regionName << endl;
238     domainDecomposition mesh
239     (
240         IOobject
241         (
242             regionName,
243             runTime.timeName(),
244             runTime
245         )
246     );
248     // Decompose the mesh
249     if (!decomposeFieldsOnly)
250     {
251         mesh.decomposeMesh(filterPatches);
253         mesh.writeDecomposition();
255         if (writeCellDist)
256         {
257             const labelList& procIds = mesh.cellToProc();
259             // Write the decomposition as labelList for use with 'manual'
260             // decomposition method.
261             labelIOList cellDecomposition
262             (
263                 IOobject
264                 (
265                     "cellDecomposition",
266                     mesh.facesInstance(),
267                     mesh,
268                     IOobject::NO_READ,
269                     IOobject::NO_WRITE,
270                     false
271                 ),
272                 procIds
273             );
274             cellDecomposition.write();
276             Info<< nl << "Wrote decomposition to "
277                 << cellDecomposition.objectPath()
278                 << " for use in manual decomposition." << endl;
280             // Write as volScalarField for post-processing
281             volScalarField cellDist
282             (
283                 IOobject
284                 (
285                     "cellDist",
286                     runTime.timeName(),
287                     mesh,
288                     IOobject::NO_READ,
289                     IOobject::NO_WRITE
290                 ),
291                 mesh,
292                 dimensionedScalar("cellDist", dimless, 0),
293                 zeroGradientFvPatchScalarField::typeName
294             );
296             forAll(procIds, celli)
297             {
298                cellDist[celli] = procIds[celli];
299             }
301             cellDist.write();
303             Info<< nl << "Wrote decomposition as volScalarField to "
304                 << cellDist.name() << " for use in post-processing."
305                 << endl;
306         }
307     }
310     // Search for list of objects for this time
311     IOobjectList objects(mesh, runTime.timeName());
313     // Construct the vol fields
314     // ~~~~~~~~~~~~~~~~~~~~~~~~
315     PtrList<volScalarField> volScalarFields;
316     readFields(mesh, objects, volScalarFields);
318     PtrList<volVectorField> volVectorFields;
319     readFields(mesh, objects, volVectorFields);
321     PtrList<volSphericalTensorField> volSphericalTensorFields;
322     readFields(mesh, objects, volSphericalTensorFields);
324     PtrList<volSymmTensorField> volSymmTensorFields;
325     readFields(mesh, objects, volSymmTensorFields);
327     PtrList<volTensorField> volTensorFields;
328     readFields(mesh, objects, volTensorFields);
331     // Construct the surface fields
332     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
333     PtrList<surfaceScalarField> surfaceScalarFields;
334     readFields(mesh, objects, surfaceScalarFields);
335     PtrList<surfaceVectorField> surfaceVectorFields;
336     readFields(mesh, objects, surfaceVectorFields);
337     PtrList<surfaceSphericalTensorField> surfaceSphericalTensorFields;
338     readFields(mesh, objects, surfaceSphericalTensorFields);
339     PtrList<surfaceSymmTensorField> surfaceSymmTensorFields;
340     readFields(mesh, objects, surfaceSymmTensorFields);
341     PtrList<surfaceTensorField> surfaceTensorFields;
342     readFields(mesh, objects, surfaceTensorFields);
345     // Construct the point fields
346     // ~~~~~~~~~~~~~~~~~~~~~~~~~~
347     pointMesh pMesh(mesh);
349     PtrList<pointScalarField> pointScalarFields;
350     readFields(pMesh, objects, pointScalarFields);
352     PtrList<pointVectorField> pointVectorFields;
353     readFields(pMesh, objects, pointVectorFields);
355     PtrList<pointSphericalTensorField> pointSphericalTensorFields;
356     readFields(pMesh, objects, pointSphericalTensorFields);
358     PtrList<pointSymmTensorField> pointSymmTensorFields;
359     readFields(pMesh, objects, pointSymmTensorFields);
361     PtrList<pointTensorField> pointTensorFields;
362     readFields(pMesh, objects, pointTensorFields);
365     // Construct the tetPoint fields
366     // ~~~~~~~~~~~~~~~~~~~~~~~~~~
367     tetPolyMesh* tetMeshPtr = NULL;
369     PtrList<tetPointScalarField> tetPointScalarFields;
370     PtrList<tetPointVectorField> tetPointVectorFields;
371     PtrList<tetPointSphericalTensorField> tetPointSphericalTensorFields;
372     PtrList<tetPointSymmTensorField> tetPointSymmTensorFields;
373     PtrList<tetPointTensorField> tetPointTensorFields;
375     PtrList<elementScalarField> elementScalarFields;
376     PtrList<elementVectorField> elementVectorFields;
378     if
379     (
380         objects.lookupClass("tetPointScalarField").size() > 0
381      || objects.lookupClass("tetPointVectorField").size() > 0
382      || objects.lookupClass("tetPointSphericalTensorField").size() > 0
383      || objects.lookupClass("tetPointSymmTensorField").size() > 0
384      || objects.lookupClass("tetPointTensorField").size() > 0
386      || objects.lookupClass("elementScalarField").size() > 0
387      || objects.lookupClass("elementVectorField").size() > 0
388     )
389     {
390         tetMeshPtr = new tetPolyMesh(mesh);
391         tetPolyMesh& tetMesh = *tetMeshPtr;
393         readFields(tetMesh, objects, tetPointScalarFields);
394         readFields(tetMesh, objects, tetPointVectorFields);
395         readFields(tetMesh, objects, tetPointSphericalTensorFields);
396         readFields(tetMesh, objects, tetPointSymmTensorFields);
397         readFields(tetMesh, objects, tetPointTensorFields);
399         readFields(tetMesh, objects, elementScalarFields);
400         readFields(tetMesh, objects, elementVectorFields);
401     }
404     // Construct the Lagrangian fields
405     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
407     fileNameList cloudDirs
408     (
409         readDir(runTime.timePath()/cloud::prefix, fileName::DIRECTORY)
410     );
412     // Particles
413     PtrList<Cloud<indexedParticle> > lagrangianPositions(cloudDirs.size());
414     // Particles per cell
415     PtrList< List<SLList<indexedParticle*>*> > cellParticles(cloudDirs.size());
417     PtrList<PtrList<labelIOField> > lagrangianLabelFields(cloudDirs.size());
418     PtrList<PtrList<scalarIOField> > lagrangianScalarFields(cloudDirs.size());
419     PtrList<PtrList<vectorIOField> > lagrangianVectorFields(cloudDirs.size());
420     PtrList<PtrList<sphericalTensorIOField> > lagrangianSphericalTensorFields
421     (
422         cloudDirs.size()
423     );
424     PtrList<PtrList<symmTensorIOField> > lagrangianSymmTensorFields
425     (
426         cloudDirs.size()
427     );
428     PtrList<PtrList<tensorIOField> > lagrangianTensorFields
429     (
430         cloudDirs.size()
431     );
433     label cloudI = 0;
435     forAll(cloudDirs, i)
436     {
437         IOobjectList sprayObjs
438         (
439             mesh,
440             runTime.timeName(),
441             cloud::prefix/cloudDirs[i]
442         );
444         IOobject* positionsPtr = sprayObjs.lookup("positions");
446         if (positionsPtr)
447         {
448             // Read lagrangian particles
449             // ~~~~~~~~~~~~~~~~~~~~~~~~~
451             Info<< "Identified lagrangian data set: " << cloudDirs[i] << endl;
453             lagrangianPositions.set
454             (
455                 cloudI,
456                 new Cloud<indexedParticle>
457                 (
458                     mesh,
459                     cloudDirs[i],
460                     false
461                 )
462             );
465             // Sort particles per cell
466             // ~~~~~~~~~~~~~~~~~~~~~~~
468             cellParticles.set
469             (
470                 cloudI,
471                 new List<SLList<indexedParticle*>*>
472                 (
473                     mesh.nCells(),
474                     static_cast<SLList<indexedParticle*>*>(NULL)
475                 )
476             );
478             label i = 0;
480             forAllIter
481             (
482                 Cloud<indexedParticle>,
483                 lagrangianPositions[cloudI],
484                 iter
485             )
486             {
487                 iter().index() = i++;
489                 label celli = iter().cell();
491                 // Check
492                 if (celli < 0 || celli >= mesh.nCells())
493                 {
494                     FatalErrorIn(args.executable())
495                         << "Illegal cell number " << celli
496                         << " for particle with index " << iter().index()
497                         << " at position " << iter().position() << nl
498                         << "Cell number should be between 0 and "
499                         << mesh.nCells()-1 << nl
500                         << "On this mesh the particle should be in cell "
501                         << mesh.findCell(iter().position())
502                         << exit(FatalError);
503                 }
505                 if (!cellParticles[cloudI][celli])
506                 {
507                     cellParticles[cloudI][celli] =
508                         new SLList<indexedParticle*>();
509                 }
511                 cellParticles[cloudI][celli]->append(&iter());
512             }
514             // Read fields
515             // ~~~~~~~~~~~
517             IOobjectList lagrangianObjects
518             (
519                 mesh,
520                 runTime.timeName(),
521                 cloud::prefix/cloudDirs[cloudI]
522             );
524             lagrangianFieldDecomposer::readFields
525             (
526                 cloudI,
527                 lagrangianObjects,
528                 lagrangianLabelFields
529             );
531             lagrangianFieldDecomposer::readFields
532             (
533                 cloudI,
534                 lagrangianObjects,
535                 lagrangianScalarFields
536             );
538             lagrangianFieldDecomposer::readFields
539             (
540                 cloudI,
541                 lagrangianObjects,
542                 lagrangianVectorFields
543             );
545             lagrangianFieldDecomposer::readFields
546             (
547                 cloudI,
548                 lagrangianObjects,
549                 lagrangianSphericalTensorFields
550             );
552             lagrangianFieldDecomposer::readFields
553             (
554                 cloudI,
555                 lagrangianObjects,
556                 lagrangianSymmTensorFields
557             );
559             lagrangianFieldDecomposer::readFields
560             (
561                 cloudI,
562                 lagrangianObjects,
563                 lagrangianTensorFields
564             );
566             cloudI++;
567         }
568     }
570     lagrangianPositions.setSize(cloudI);
571     cellParticles.setSize(cloudI);
572     lagrangianLabelFields.setSize(cloudI);
573     lagrangianScalarFields.setSize(cloudI);
574     lagrangianVectorFields.setSize(cloudI);
575     lagrangianSphericalTensorFields.setSize(cloudI);
576     lagrangianSymmTensorFields.setSize(cloudI);
577     lagrangianTensorFields.setSize(cloudI);
580     // Any uniform data to copy/link?
581     fileName uniformDir("uniform");
583     if (isDir(runTime.timePath()/uniformDir))
584     {
585         Info<< "Detected additional non-decomposed files in "
586             << runTime.timePath()/uniformDir
587             << endl;
588     }
589     else
590     {
591         uniformDir.clear();
592     }
594     Info<< endl;
596     // Split the fields over processors
597     for (label procI = 0; procI < mesh.nProcs(); procI++)
598     {
599         Info<< "Processor " << procI << ": field transfer" << endl;
601         // open the database
602         Time processorDb
603         (
604             Time::controlDictName,
605             args.rootPath(),
606             args.caseName()/fileName(word("processor") + name(procI))
607         );
609         processorDb.setTime(runTime);
611         // Remove files remnants that can cause horrible problems
612         // - mut and nut are used to mark the new turbulence models,
613         //   their existence prevents old models from being upgraded
614         // 1.6.x merge.  HJ, 25/Aug/2010
615         {
616             fileName timeDir(processorDb.path()/processorDb.timeName());
618             rm(timeDir/"mut");
619             rm(timeDir/"nut");
620             rm(timeDir/"mut.gz");
621             rm(timeDir/"nut.gz");
622         }
624         // read the mesh
625         fvMesh procMesh
626         (
627             IOobject
628             (
629                 regionName,
630                 processorDb.timeName(),
631                 processorDb
632             )
633         );
635         labelIOList cellProcAddressing
636         (
637             IOobject
638             (
639                 "cellProcAddressing",
640                 procMesh.facesInstance(),
641                 procMesh.meshSubDir,
642                 procMesh,
643                 IOobject::MUST_READ,
644                 IOobject::NO_WRITE
645             )
646         );
648         labelIOList boundaryProcAddressing
649         (
650             IOobject
651             (
652                 "boundaryProcAddressing",
653                 procMesh.facesInstance(),
654                 procMesh.meshSubDir,
655                 procMesh,
656                 IOobject::MUST_READ,
657                 IOobject::NO_WRITE
658             )
659         );
661         // FV fields
662         if
663         (
664             volScalarFields.size()
665          || volVectorFields.size()
666          || volSphericalTensorFields.size()
667          || volSymmTensorFields.size()
668          || volTensorFields.size()
669          || surfaceScalarFields.size()
670          || surfaceVectorFields.size()
671          || surfaceSphericalTensorFields.size()
672          || surfaceSymmTensorFields.size()
673          || surfaceTensorFields.size()
674         )
675         {
676             labelIOList faceProcAddressing
677             (
678                 IOobject
679                 (
680                     "faceProcAddressing",
681                     procMesh.facesInstance(),
682                     procMesh.meshSubDir,
683                     procMesh,
684                     IOobject::MUST_READ,
685                     IOobject::NO_WRITE
686                 )
687             );
689             fvFieldDecomposer fieldDecomposer
690             (
691                 mesh,
692                 procMesh,
693                 faceProcAddressing,
694                 cellProcAddressing,
695                 boundaryProcAddressing
696             );
698             fieldDecomposer.decomposeFields(volScalarFields);
699             fieldDecomposer.decomposeFields(volVectorFields);
700             fieldDecomposer.decomposeFields(volSphericalTensorFields);
701             fieldDecomposer.decomposeFields(volSymmTensorFields);
702             fieldDecomposer.decomposeFields(volTensorFields);
704             fieldDecomposer.decomposeFields(surfaceScalarFields);
705             fieldDecomposer.decomposeFields(surfaceVectorFields);
706             fieldDecomposer.decomposeFields(surfaceSphericalTensorFields);
707             fieldDecomposer.decomposeFields(surfaceSymmTensorFields);
708             fieldDecomposer.decomposeFields(surfaceTensorFields);
709         }
712         // Point fields
713         if
714         (
715             pointScalarFields.size()
716          || pointVectorFields.size()
717          || pointSphericalTensorFields.size()
718          || pointSymmTensorFields.size()
719          || pointTensorFields.size()
720         )
721         {
722             labelIOList pointProcAddressing
723             (
724                 IOobject
725                 (
726                     "pointProcAddressing",
727                     procMesh.facesInstance(),
728                     procMesh.meshSubDir,
729                     procMesh,
730                     IOobject::MUST_READ,
731                     IOobject::NO_WRITE
732                 )
733             );
735             pointMesh procPMesh(procMesh, true);
737             pointFieldDecomposer fieldDecomposer
738             (
739                 pMesh,
740                 procPMesh,
741                 pointProcAddressing,
742                 boundaryProcAddressing
743             );
745             fieldDecomposer.decomposeFields(pointScalarFields);
746             fieldDecomposer.decomposeFields(pointVectorFields);
747             fieldDecomposer.decomposeFields(pointSphericalTensorFields);
748             fieldDecomposer.decomposeFields(pointSymmTensorFields);
749             fieldDecomposer.decomposeFields(pointTensorFields);
750         }
753         // tetPoint fields
754         if (tetMeshPtr)
755         {
756             const tetPolyMesh& tetMesh = *tetMeshPtr;
757             tetPolyMesh procTetMesh(procMesh);
759             // Read the point addressing information
760             labelIOList pointProcAddressing
761             (
762                 IOobject
763                 (
764                     "pointProcAddressing",
765                     procMesh.facesInstance(),
766                     procMesh.meshSubDir,
767                     procMesh,
768                     IOobject::MUST_READ,
769                     IOobject::NO_WRITE
770                 )
771             );
773             // Read the point addressing information
774             labelIOList faceProcAddressing
775             (
776                 IOobject
777                 (
778                     "faceProcAddressing",
779                     procMesh.facesInstance(),
780                     procMesh.meshSubDir,
781                     procMesh,
782                     IOobject::MUST_READ,
783                     IOobject::NO_WRITE
784                 )
785             );
787             tetPointFieldDecomposer fieldDecomposer
788             (
789                 tetMesh,
790                 procTetMesh,
791                 pointProcAddressing,
792                 faceProcAddressing,
793                 cellProcAddressing,
794                 boundaryProcAddressing
795             );
797             fieldDecomposer.decomposeFields(tetPointScalarFields);
798             fieldDecomposer.decomposeFields(tetPointVectorFields);
799             fieldDecomposer.decomposeFields(tetPointSphericalTensorFields);
800             fieldDecomposer.decomposeFields(tetPointSymmTensorFields);
801             fieldDecomposer.decomposeFields(tetPointTensorFields);
803             fieldDecomposer.decomposeFields(elementScalarFields);
804             fieldDecomposer.decomposeFields(elementVectorFields);
805         }
808         // If there is lagrangian data write it out
809         forAll(lagrangianPositions, cloudI)
810         {
811             if (lagrangianPositions[cloudI].size())
812             {
813                 lagrangianFieldDecomposer fieldDecomposer
814                 (
815                     mesh,
816                     procMesh,
817                     cellProcAddressing,
818                     cloudDirs[cloudI],
819                     lagrangianPositions[cloudI],
820                     cellParticles[cloudI]
821                 );
823                 // Lagrangian fields
824                 if
825                 (
826                     lagrangianLabelFields[cloudI].size()
827                  || lagrangianScalarFields[cloudI].size()
828                  || lagrangianVectorFields[cloudI].size()
829                  || lagrangianSphericalTensorFields[cloudI].size()
830                  || lagrangianSymmTensorFields[cloudI].size()
831                  || lagrangianTensorFields[cloudI].size()
832                 )
833                 {
834                     fieldDecomposer.decomposeFields
835                     (
836                         cloudDirs[cloudI],
837                         lagrangianLabelFields[cloudI]
838                     );
839                     fieldDecomposer.decomposeFields
840                     (
841                         cloudDirs[cloudI],
842                         lagrangianScalarFields[cloudI]
843                     );
844                     fieldDecomposer.decomposeFields
845                     (
846                         cloudDirs[cloudI],
847                         lagrangianVectorFields[cloudI]
848                     );
849                     fieldDecomposer.decomposeFields
850                     (
851                         cloudDirs[cloudI],
852                         lagrangianSphericalTensorFields[cloudI]
853                     );
854                     fieldDecomposer.decomposeFields
855                     (
856                         cloudDirs[cloudI],
857                         lagrangianSymmTensorFields[cloudI]
858                     );
859                     fieldDecomposer.decomposeFields
860                     (
861                         cloudDirs[cloudI],
862                         lagrangianTensorFields[cloudI]
863                     );
864                 }
865             }
866         }
869         // Any non-decomposed data to copy?
870         if (uniformDir.size())
871         {
872             const fileName timePath = processorDb.timePath();
874             if (copyUniform || mesh.distributed())
875             {
876                 cp
877                 (
878                     runTime.timePath()/uniformDir,
879                     timePath/uniformDir
880                 );
881             }
882             else
883             {
884                 // Link with relative paths
885                 const string parentPath = string("..")/"..";
887                 fileName currentDir(cwd());
888                 chDir(timePath);
889                 if (!exists(uniformDir))
890                 {
891                     ln
892                     (
893                         parentPath/runTime.timeName()/uniformDir,
894                         uniformDir
895                     );
896                     chDir(currentDir);
897                 }
898             }
899         }
900     }
903     if (tetMeshPtr)
904     {
905         delete tetMeshPtr;
906         tetMeshPtr = NULL;
907     }
910     // Finite area mesh and field decomposition
912     IOobject faMeshBoundaryIOobj
913     (
914         "faBoundary",
915         mesh.time().findInstance
916         (
917             mesh.dbDir()/polyMesh::meshSubDir,
918             "boundary"
919         ),
920         faMesh::meshSubDir,
921         mesh,
922         IOobject::READ_IF_PRESENT,
923         IOobject::NO_WRITE
924     );
927     if(faMeshBoundaryIOobj.headerOk())
928     {
929         Info << "\nFinite area mesh decomposition" << endl;
931         faMeshDecomposition aMesh(mesh);
933         aMesh.decomposeMesh(filterPatches);
935         aMesh.writeDecomposition();
938         // Construct the area fields
939         // ~~~~~~~~~~~~~~~~~~~~~~~~
940         PtrList<areaScalarField> areaScalarFields;
941         readFields(aMesh, objects, areaScalarFields);
943         PtrList<areaVectorField> areaVectorFields;
944         readFields(aMesh, objects, areaVectorFields);
946         PtrList<areaSphericalTensorField> areaSphericalTensorFields;
947         readFields(aMesh, objects, areaSphericalTensorFields);
949         PtrList<areaSymmTensorField> areaSymmTensorFields;
950         readFields(aMesh, objects, areaSymmTensorFields);
952         PtrList<areaTensorField> areaTensorFields;
953         readFields(aMesh, objects, areaTensorFields);
956         // Construct the edge fields
957         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
958         PtrList<edgeScalarField> edgeScalarFields;
959         readFields(aMesh, objects, edgeScalarFields);
961         Info << endl;
963         // Split the fields over processors
964         for (label procI = 0; procI < mesh.nProcs(); procI++)
965         {
966             Info<< "Processor " << procI
967                 << ": finite area field transfer" << endl;
969             // open the database
970             Time processorDb
971             (
972                 Time::controlDictName,
973                 args.rootPath(),
974                 args.caseName()/fileName(word("processor") + name(procI))
975             );
977             processorDb.setTime(runTime);
979             // Read the mesh
980             fvMesh procFvMesh
981             (
982                 IOobject
983                 (
984                     regionName,
985                     processorDb.timeName(),
986                     processorDb
987                 )
988             );
990             faMesh procMesh(procFvMesh);
992             labelIOList faceProcAddressing
993             (
994                 IOobject
995                 (
996                     "faceProcAddressing",
997                     "constant",
998                     procMesh.meshSubDir,
999                     procFvMesh,
1000                     IOobject::MUST_READ,
1001                     IOobject::NO_WRITE
1002                 )
1003             );
1005             labelIOList boundaryProcAddressing
1006             (
1007                 IOobject
1008                 (
1009                     "boundaryProcAddressing",
1010                     "constant",
1011                     procMesh.meshSubDir,
1012                     procFvMesh,
1013                     IOobject::MUST_READ,
1014                     IOobject::NO_WRITE
1015                 )
1016             );
1018             // FA fields
1019             if
1020             (
1021                 areaScalarFields.size()
1022              || areaVectorFields.size()
1023              || areaSphericalTensorFields.size()
1024              || areaSymmTensorFields.size()
1025              || areaTensorFields.size()
1026              || edgeScalarFields.size()
1027             )
1028             {
1029                 labelIOList edgeProcAddressing
1030                 (
1031                     IOobject
1032                     (
1033                         "edgeProcAddressing",
1034                         "constant",
1035                         procMesh.meshSubDir,
1036                         procFvMesh,
1037                         IOobject::MUST_READ,
1038                         IOobject::NO_WRITE
1039                     )
1040                 );
1042                 faFieldDecomposer fieldDecomposer
1043                 (
1044                     aMesh,
1045                     procMesh,
1046                     edgeProcAddressing,
1047                     faceProcAddressing,
1048                     boundaryProcAddressing
1049                 );
1051                 fieldDecomposer.decomposeFields(areaScalarFields);
1052                 fieldDecomposer.decomposeFields(areaVectorFields);
1053                 fieldDecomposer.decomposeFields(areaSphericalTensorFields);
1054                 fieldDecomposer.decomposeFields(areaSymmTensorFields);
1055                 fieldDecomposer.decomposeFields(areaTensorFields);
1057                 fieldDecomposer.decomposeFields(edgeScalarFields);
1058             }
1059         }
1060     }
1063     Info<< "\nEnd.\n" << endl;
1065     return 0;
1069 // ************************************************************************* //