00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070 #include <OpenFOAM/argList.H>
00071 #include <OpenFOAM/timeSelector.H>
00072
00073 #include <finiteVolume/fvCFD.H>
00074 #include <OpenFOAM/IOobjectList.H>
00075 #include "processorMeshes.H"
00076 #include "fvFieldReconstructor.H"
00077 #include "pointFieldReconstructor.H"
00078 #include "reconstructLagrangian.H"
00079
00080
00081
00082 int main(int argc, char *argv[])
00083 {
00084
00085
00086 timeSelector::addOptions(true, true);
00087 argList::noParallel();
00088 # include <OpenFOAM/addRegionOption.H>
00089 argList::validOptions.insert("fields", "\"(list of fields)\"");
00090 argList::validOptions.insert("noLagrangian", "");
00091
00092 # include <OpenFOAM/setRootCase.H>
00093 # include <OpenFOAM/createTime.H>
00094
00095 HashSet<word> selectedFields;
00096 if (args.optionFound("fields"))
00097 {
00098 args.optionLookup("fields")() >> selectedFields;
00099 }
00100
00101 bool noLagrangian = args.optionFound("noLagrangian");
00102
00103
00104 label nProcs = 0;
00105 while (isDir(args.path()/(word("processor") + name(nProcs))))
00106 {
00107 ++nProcs;
00108 }
00109
00110 if (!nProcs)
00111 {
00112 FatalErrorIn(args.executable())
00113 << "No processor* directories found"
00114 << exit(FatalError);
00115 }
00116
00117
00118 PtrList<Time> databases(nProcs);
00119
00120 forAll (databases, procI)
00121 {
00122 databases.set
00123 (
00124 procI,
00125 new Time
00126 (
00127 Time::controlDictName,
00128 args.rootPath(),
00129 args.caseName()/fileName(word("processor") + name(procI))
00130 )
00131 );
00132 }
00133
00134
00135
00136 instantList timeDirs = timeSelector::select
00137 (
00138 databases[0].times(),
00139 args
00140 );
00141
00142 if (timeDirs.empty())
00143 {
00144 FatalErrorIn(args.executable())
00145 << "No times selected"
00146 << exit(FatalError);
00147 }
00148
00149 # include <OpenFOAM/createNamedMesh.H>
00150 fileName regionPrefix = "";
00151 if (regionName != fvMesh::defaultRegion)
00152 {
00153 regionPrefix = regionName;
00154 }
00155
00156
00157 forAll (databases, procI)
00158 {
00159 databases[procI].setTime(runTime.timeName(), runTime.timeIndex());
00160 }
00161
00162
00163 processorMeshes procMeshes(databases, regionName);
00164
00165
00166
00167
00168 # include "checkFaceAddressingComp.H"
00169
00170
00171 forAll (timeDirs, timeI)
00172 {
00173
00174 runTime.setTime(timeDirs[timeI], timeI);
00175
00176 Info << "Time = " << runTime.timeName() << endl << endl;
00177
00178
00179 forAll (databases, procI)
00180 {
00181 databases[procI].setTime(timeDirs[timeI], timeI);
00182 }
00183
00184
00185 fvMesh::readUpdateState meshStat = mesh.readUpdate();
00186
00187 fvMesh::readUpdateState procStat = procMeshes.readUpdate();
00188
00189 if (procStat == fvMesh::POINTS_MOVED)
00190 {
00191
00192 procMeshes.reconstructPoints(mesh);
00193 }
00194 else if (meshStat != procStat)
00195 {
00196 WarningIn(args.executable())
00197 << "readUpdate for the reconstructed mesh:" << meshStat << nl
00198 << "readUpdate for the processor meshes :" << procStat << nl
00199 << "These should be equal or your addressing"
00200 << " might be incorrect."
00201 << " Please check your time directories for any "
00202 << "mesh directories." << endl;
00203 }
00204
00205
00206
00207 IOobjectList objects(procMeshes.meshes()[0], databases[0].timeName());
00208
00209
00210
00211
00212 if
00213 (
00214 objects.lookupClass(volScalarField::typeName).size()
00215 || objects.lookupClass(volVectorField::typeName).size()
00216 || objects.lookupClass(volSphericalTensorField::typeName).size()
00217 || objects.lookupClass(volSymmTensorField::typeName).size()
00218 || objects.lookupClass(volTensorField::typeName).size()
00219 || objects.lookupClass(surfaceScalarField::typeName).size()
00220 || objects.lookupClass(surfaceVectorField::typeName).size()
00221 || objects.lookupClass(surfaceSphericalTensorField::typeName).size()
00222 || objects.lookupClass(surfaceSymmTensorField::typeName).size()
00223 || objects.lookupClass(surfaceTensorField::typeName).size()
00224 )
00225 {
00226 Info << "Reconstructing FV fields" << nl << endl;
00227
00228 fvFieldReconstructor fvReconstructor
00229 (
00230 mesh,
00231 procMeshes.meshes(),
00232 procMeshes.faceProcAddressing(),
00233 procMeshes.cellProcAddressing(),
00234 procMeshes.boundaryProcAddressing()
00235 );
00236
00237 fvReconstructor.reconstructFvVolumeFields<scalar>
00238 (
00239 objects,
00240 selectedFields
00241 );
00242 fvReconstructor.reconstructFvVolumeFields<vector>
00243 (
00244 objects,
00245 selectedFields
00246 );
00247 fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
00248 (
00249 objects,
00250 selectedFields
00251 );
00252 fvReconstructor.reconstructFvVolumeFields<symmTensor>
00253 (
00254 objects,
00255 selectedFields
00256 );
00257 fvReconstructor.reconstructFvVolumeFields<tensor>
00258 (
00259 objects,
00260 selectedFields
00261 );
00262
00263 fvReconstructor.reconstructFvSurfaceFields<scalar>
00264 (
00265 objects,
00266 selectedFields
00267 );
00268 fvReconstructor.reconstructFvSurfaceFields<vector>
00269 (
00270 objects,
00271 selectedFields
00272 );
00273 fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
00274 (
00275 objects,
00276 selectedFields
00277 );
00278 fvReconstructor.reconstructFvSurfaceFields<symmTensor>
00279 (
00280 objects,
00281 selectedFields
00282 );
00283 fvReconstructor.reconstructFvSurfaceFields<tensor>
00284 (
00285 objects,
00286 selectedFields
00287 );
00288 }
00289 else
00290 {
00291 Info << "No FV fields" << nl << endl;
00292 }
00293
00294
00295
00296 if
00297 (
00298 objects.lookupClass(pointScalarField::typeName).size()
00299 || objects.lookupClass(pointVectorField::typeName).size()
00300 || objects.lookupClass(pointSphericalTensorField::typeName).size()
00301 || objects.lookupClass(pointSymmTensorField::typeName).size()
00302 || objects.lookupClass(pointTensorField::typeName).size()
00303 )
00304 {
00305 Info << "Reconstructing point fields" << nl << endl;
00306
00307 pointMesh pMesh(mesh);
00308 PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
00309
00310 forAll (pMeshes, procI)
00311 {
00312 pMeshes.set(procI, new pointMesh(procMeshes.meshes()[procI]));
00313 }
00314
00315 pointFieldReconstructor pointReconstructor
00316 (
00317 pMesh,
00318 pMeshes,
00319 procMeshes.pointProcAddressing(),
00320 procMeshes.boundaryProcAddressing()
00321 );
00322
00323 pointReconstructor.reconstructFields<scalar>(objects);
00324 pointReconstructor.reconstructFields<vector>(objects);
00325 pointReconstructor.reconstructFields<sphericalTensor>(objects);
00326 pointReconstructor.reconstructFields<symmTensor>(objects);
00327 pointReconstructor.reconstructFields<tensor>(objects);
00328 }
00329 else
00330 {
00331 Info << "No point fields" << nl << endl;
00332 }
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342 if (!noLagrangian)
00343 {
00344 HashTable<IOobjectList> cloudObjects;
00345
00346 forAll (databases, procI)
00347 {
00348 fileNameList cloudDirs
00349 (
00350 readDir
00351 (
00352 databases[procI].timePath()/regionPrefix/cloud::prefix,
00353 fileName::DIRECTORY
00354 )
00355 );
00356
00357 forAll (cloudDirs, i)
00358 {
00359
00360 HashTable<IOobjectList>::const_iterator iter =
00361 cloudObjects.find(cloudDirs[i]);
00362
00363 if (iter == cloudObjects.end())
00364 {
00365
00366 IOobjectList sprayObjs
00367 (
00368 procMeshes.meshes()[procI],
00369 databases[procI].timeName(),
00370 cloud::prefix/cloudDirs[i]
00371 );
00372
00373 IOobject* positionsPtr = sprayObjs.lookup("positions");
00374
00375 if (positionsPtr)
00376 {
00377 cloudObjects.insert(cloudDirs[i], sprayObjs);
00378 }
00379 }
00380 }
00381 }
00382
00383
00384 if (cloudObjects.size())
00385 {
00386
00387 forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
00388 {
00389 const word cloudName = string::validate<word>(iter.key());
00390
00391
00392 const IOobjectList& sprayObjs = iter();
00393
00394 Info<< "Reconstructing lagrangian fields for cloud "
00395 << cloudName << nl << endl;
00396
00397 reconstructLagrangianPositions
00398 (
00399 mesh,
00400 cloudName,
00401 procMeshes.meshes(),
00402 procMeshes.faceProcAddressing(),
00403 procMeshes.cellProcAddressing()
00404 );
00405 reconstructLagrangianFields<label>
00406 (
00407 cloudName,
00408 mesh,
00409 procMeshes.meshes(),
00410 sprayObjs
00411 );
00412 reconstructLagrangianFields<scalar>
00413 (
00414 cloudName,
00415 mesh,
00416 procMeshes.meshes(),
00417 sprayObjs
00418 );
00419 reconstructLagrangianFields<vector>
00420 (
00421 cloudName,
00422 mesh,
00423 procMeshes.meshes(),
00424 sprayObjs
00425 );
00426 reconstructLagrangianFields<sphericalTensor>
00427 (
00428 cloudName,
00429 mesh,
00430 procMeshes.meshes(),
00431 sprayObjs
00432 );
00433 reconstructLagrangianFields<symmTensor>
00434 (
00435 cloudName,
00436 mesh,
00437 procMeshes.meshes(),
00438 sprayObjs
00439 );
00440 reconstructLagrangianFields<tensor>
00441 (
00442 cloudName,
00443 mesh,
00444 procMeshes.meshes(),
00445 sprayObjs
00446 );
00447 }
00448 }
00449 else
00450 {
00451 Info << "No lagrangian fields" << nl << endl;
00452 }
00453 }
00454
00455
00456
00457
00458 fileName uniformDir0 = databases[0].timePath()/"uniform";
00459 if (isDir(uniformDir0))
00460 {
00461 cp(uniformDir0, runTime.timePath());
00462 }
00463 }
00464
00465 Info<< "End.\n" << endl;
00466
00467 return 0;
00468 }
00469
00470
00471