Go to the documentation of this file.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 #include "MapLagrangianFields.H"
00027 #include <lagrangian/Cloud.H>
00028 #include <lagrangian/passiveParticle.H>
00029 #include <meshTools/meshSearch.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035
00036 static const scalar perturbFactor = 1E-6;
00037
00038
00039
00040
00041 static label findCell(const meshSearch& meshSearcher, const point& pt)
00042 {
00043 const polyMesh& mesh = meshSearcher.mesh();
00044
00045
00046 label cellI = meshSearcher.findCell(pt);
00047
00048 if (cellI >= 0)
00049 {
00050 return cellI;
00051 }
00052 else
00053 {
00054
00055
00056
00057 label faceI = meshSearcher.findNearestBoundaryFace(pt);
00058
00059 if (faceI >= 0)
00060 {
00061 const point& cc = mesh.cellCentres()[mesh.faceOwner()[faceI]];
00062 const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
00063
00064 return meshSearcher.findCell(perturbPt);
00065 }
00066 }
00067 return -1;
00068 }
00069
00070
00071 void mapLagrangian(const meshToMesh& meshToMeshInterp)
00072 {
00073
00074
00075
00076
00077 const labelList& cellAddressing = meshToMeshInterp.cellAddressing();
00078
00079
00080
00081
00082 labelListList sourceToTargets
00083 (
00084 invertOneToMany(meshToMeshInterp.fromMesh().nCells(), cellAddressing)
00085 );
00086
00087 const fvMesh& meshSource = meshToMeshInterp.fromMesh();
00088 const fvMesh& meshTarget = meshToMeshInterp.toMesh();
00089 const pointField& targetCc = meshTarget.cellCentres();
00090
00091
00092 fileNameList cloudDirs
00093 (
00094 readDir
00095 (
00096 meshSource.time().timePath()/cloud::prefix,
00097 fileName::DIRECTORY
00098 )
00099 );
00100
00101 forAll(cloudDirs, cloudI)
00102 {
00103
00104 IOobjectList objects
00105 (
00106 meshSource,
00107 meshSource.time().timeName(),
00108 cloud::prefix/cloudDirs[cloudI]
00109 );
00110
00111 IOobject* positionsPtr = objects.lookup("positions");
00112
00113 if (positionsPtr)
00114 {
00115 Info<< nl << " processing cloud " << cloudDirs[cloudI] << endl;
00116
00117
00118 Cloud<passiveParticle> sourceParcels
00119 (
00120 meshSource,
00121 cloudDirs[cloudI],
00122 false
00123 );
00124 Info<< " read " << sourceParcels.size()
00125 << " parcels from source mesh." << endl;
00126
00127
00128 Cloud<passiveParticle> targetParcels
00129 (
00130 meshTarget,
00131 cloudDirs[cloudI],
00132 IDLList<passiveParticle>()
00133 );
00134
00135 label sourceParticleI = 0;
00136
00137
00138 DynamicList<label> addParticles(sourceParcels.size());
00139
00140
00141 labelHashSet unmappedSource(sourceParcels.size());
00142
00143
00144
00145
00146
00147
00148
00149 forAllConstIter(Cloud<passiveParticle>, sourceParcels, iter)
00150 {
00151 bool foundCell = false;
00152
00153
00154 if (iter().cell() >= 0)
00155 {
00156 const labelList& targetCells =
00157 sourceToTargets[iter().cell()];
00158
00159
00160
00161
00162
00163 forAll(targetCells, i)
00164 {
00165
00166 autoPtr<passiveParticle> newPtr
00167 (
00168 new passiveParticle
00169 (
00170 targetParcels,
00171 targetCc[targetCells[i]],
00172 targetCells[i]
00173 )
00174 );
00175 passiveParticle& newP = newPtr();
00176
00177 scalar fraction = 0;
00178 label faceI = newP.track(iter().position(), fraction);
00179
00180 if (faceI < 0 && newP.cell() >= 0)
00181 {
00182
00183 foundCell = true;
00184 addParticles.append(sourceParticleI);
00185 targetParcels.addParticle(newPtr.ptr());
00186 break;
00187 }
00188 }
00189 }
00190
00191 if (!foundCell)
00192 {
00193
00194 unmappedSource.insert(sourceParticleI);
00195 }
00196
00197 sourceParticleI++;
00198 }
00199
00200 Info<< " after meshToMesh addressing found "
00201 << targetParcels.size()
00202 << " parcels in target mesh." << endl;
00203
00204
00205
00206
00207
00208 if (unmappedSource.size())
00209 {
00210 meshSearch targetSearcher(meshTarget, false);
00211
00212 sourceParticleI = 0;
00213
00214 forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
00215 {
00216 if (unmappedSource.found(sourceParticleI))
00217 {
00218 label targetCell =
00219 findCell(targetSearcher, iter().position());
00220
00221 if (targetCell >= 0)
00222 {
00223 unmappedSource.erase(sourceParticleI);
00224 addParticles.append(sourceParticleI);
00225 iter().cell()=targetCell;
00226 targetParcels.addParticle
00227 (
00228 sourceParcels.remove(&iter())
00229 );
00230 }
00231 }
00232 sourceParticleI++;
00233 }
00234 }
00235 addParticles.shrink();
00236
00237 Info<< " after additional mesh searching found "
00238 << targetParcels.size() << " parcels in target mesh." << endl;
00239
00240 if (addParticles.size())
00241 {
00242 IOPosition<passiveParticle>(targetParcels).write();
00243
00244
00245
00246
00247
00248
00249
00250 MapLagrangianFields<label>
00251 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00252 MapLagrangianFields<scalar>
00253 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00254 MapLagrangianFields<vector>
00255 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00256 MapLagrangianFields<sphericalTensor>
00257 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00258 MapLagrangianFields<symmTensor>
00259 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00260 MapLagrangianFields<tensor>
00261 (cloudDirs[cloudI], objects, meshToMeshInterp, addParticles);
00262 }
00263 }
00264 }
00265 }
00266
00267
00268
00269
00270 }
00271
00272