Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
Henry
committed
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "mappedPatchBase.H"
#include "addToRunTimeSelectionTable.H"
#include "ListListOps.H"
#include "meshSearchMeshObject.H"
#include "meshTools.H"
#include "OFstream.H"
#include "Random.H"
#include "treeDataFace.H"
#include "treeDataPoint.H"
#include "indexedOctree.H"
#include "polyMesh.H"
#include "polyPatch.H"
#include "Time.H"
#include "mapDistribute.H"
#include "SubField.H"
#include "syncTools.H"
#include "treeDataCell.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(mappedPatchBase, 0);
template<>
const char* Foam::NamedEnum
<
Foam::mappedPatchBase::sampleMode,
>::names[] =
{
"nearestCell",
"nearestPatchFace",
"nearestFace",
"nearestOnlyCell"
};
template<>
const char* Foam::NamedEnum
<
Foam::mappedPatchBase::offsetMode,
3
>::names[] =
{
"uniform",
"nonuniform",
"normal"
};
}
const Foam::NamedEnum<Foam::mappedPatchBase::sampleMode, 6>
Foam::mappedPatchBase::sampleModeNames_;
const Foam::NamedEnum<Foam::mappedPatchBase::offsetMode, 3>
Foam::mappedPatchBase::offsetModeNames_;
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::pointField> Foam::mappedPatchBase::facePoints
(
const polyPatch& pp
) const
{
const polyMesh& mesh = pp.boundaryMesh().mesh();
// Force construction of min-tet decomp
(void)mesh.tetBasePtIs();
// Initialise to face-centre
tmp<pointField> tfacePoints(new pointField(patch_.size()));
pointField& facePoints = tfacePoints();
forAll(pp, faceI)
{
facePoints[faceI] = facePoint
(
mesh,
pp.start()+faceI,
Henry
committed
polyMesh::FACE_DIAG_TRIS
}
return tfacePoints;
}
void Foam::mappedPatchBase::collectSamples
(
const pointField& facePoints,
pointField& samples,
labelList& patchFaceProcs,
labelList& patchFaces,
pointField& patchFc
) const
{
// Collect all sample points and the faces they come from.
{
List<pointField> globalFc(Pstream::nProcs());
globalFc[Pstream::myProcNo()] = facePoints;
Pstream::gatherList(globalFc);
Pstream::scatterList(globalFc);
// Rework into straight list
patchFc = ListListOps::combine<pointField>
(
globalFc,
accessOp<pointField>()
);
}
{
List<pointField> globalSamples(Pstream::nProcs());
globalSamples[Pstream::myProcNo()] = samplePoints(facePoints);
Pstream::gatherList(globalSamples);
Pstream::scatterList(globalSamples);
// Rework into straight list
samples = ListListOps::combine<pointField>
(
globalSamples,
accessOp<pointField>()
);
}
{
labelListList globalFaces(Pstream::nProcs());
globalFaces[Pstream::myProcNo()] = identity(patch_.size());
// Distribute to all processors
Pstream::gatherList(globalFaces);
Pstream::scatterList(globalFaces);
patchFaces = ListListOps::combine<labelList>
(
globalFaces,
accessOp<labelList>()
labelList nPerProc(Pstream::nProcs());
nPerProc[Pstream::myProcNo()] = patch_.size();
Pstream::gatherList(nPerProc);
Pstream::scatterList(nPerProc);
patchFaceProcs.setSize(patchFaces.size());
label sampleI = 0;
forAll(nPerProc, procI)
for (label i = 0; i < nPerProc[procI]; i++)
{
patchFaceProcs[sampleI++] = procI;
}
}
}
}
// Find the processor/cell containing the samples. Does not account
// for samples being found in two processors.
void Foam::mappedPatchBase::findSamples
(
const pointField& samples,
labelList& sampleProcs,
labelList& sampleIndices,
pointField& sampleLocations
) const
{
// Lookup the correct region
const polyMesh& mesh = sampleMesh();
// All the info for nearest. Construct to miss
List<nearInfo> nearest(samples.size());
{
case NEARESTCELL:
{
andy
committed
if (samplePatch_.size() && samplePatch_ != "none")
{
FatalErrorIn
(
"mappedPatchBase::findSamples(const pointField&,"
" labelList&, labelList&, pointField&) const"
) << "No need to supply a patch name when in "
<< sampleModeNames_[mode] << " mode." << exit(FatalError);
//- Note: face-diagonal decomposition
const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
forAll(samples, sampleI)
{
const point& sample = samples[sampleI];
label cellI = tree.findInside(sample);
if (cellI == -1)
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
else
{
const point& cc = mesh.cellCentres()[cellI];
nearest[sampleI].first() = pointIndexHit
(
true,
cc,
cellI
);
nearest[sampleI].second().first() = magSqr(cc-sample);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
}
break;
}
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
case NEARESTONLYCELL:
{
if (samplePatch_.size() && samplePatch_ != "none")
{
FatalErrorIn
(
"mappedPatchBase::findSamples(const pointField&,"
" labelList&, labelList&, pointField&) const"
) << "No need to supply a patch name when in "
<< sampleModeNames_[mode] << " mode." << exit(FatalError);
}
//- Note: face-diagonal decomposition
const indexedOctree<Foam::treeDataCell>& tree = mesh.cellTree();
forAll(samples, sampleI)
{
const point& sample = samples[sampleI];
nearest[sampleI].first() = tree.findNearest(sample, sqr(GREAT));
nearest[sampleI].second().first() = magSqr
(
nearest[sampleI].first().hitPoint()
-sample
);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
break;
}
case NEARESTPATCHFACE:
{
Random rndGen(123456);
const polyPatch& pp = samplePolyPatch();
if (pp.empty())
{
forAll(samples, sampleI)
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
}
else
{
// patch faces
const labelList patchFaces(identity(pp.size()) + pp.start());
treeBoundBox patchBb
(
treeBoundBox(pp.points(), pp.meshPoints()).extend
(
rndGen,
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
)
);
patchBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
patchBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
indexedOctree<treeDataFace> boundaryTree
(
treeDataFace // all information needed to search faces
(
false, // do not cache bb
mesh,
patchFaces // boundary faces only
),
patchBb, // overall search domain
8, // maxLevel
10, // leafsize
3.0 // duplicity
);
forAll(samples, sampleI)
{
const point& sample = samples[sampleI];
pointIndexHit& nearInfo = nearest[sampleI].first();
nearInfo = boundaryTree.findNearest
(
sample,
magSqr(patchBb.span())
);
if (!nearInfo.hit())
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() =
Pstream::myProcNo();
}
else
{
point fc(pp[nearInfo.index()].centre(pp.points()));
nearInfo.setPoint(fc);
nearest[sampleI].second().first() = magSqr(fc-sample);
nearest[sampleI].second().second() =
Pstream::myProcNo();
}
}
}
break;
}
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
case NEARESTPATCHPOINT:
{
Random rndGen(123456);
const polyPatch& pp = samplePolyPatch();
if (pp.empty())
{
forAll(samples, sampleI)
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
}
else
{
// patch (local) points
treeBoundBox patchBb
(
treeBoundBox(pp.points(), pp.meshPoints()).extend
(
rndGen,
1e-4
)
);
patchBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
patchBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
indexedOctree<treeDataPoint> boundaryTree
(
treeDataPoint // all information needed to search faces
(
mesh.points(),
pp.meshPoints() // selection of points to search on
),
patchBb, // overall search domain
8, // maxLevel
10, // leafsize
3.0 // duplicity
);
forAll(samples, sampleI)
{
const point& sample = samples[sampleI];
pointIndexHit& nearInfo = nearest[sampleI].first();
nearInfo = boundaryTree.findNearest
(
sample,
magSqr(patchBb.span())
);
if (!nearInfo.hit())
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() =
Pstream::myProcNo();
}
else
{
const point& pt = nearInfo.hitPoint();
nearest[sampleI].second().first() = magSqr(pt-sample);
nearest[sampleI].second().second() =
Pstream::myProcNo();
}
}
}
break;
}
case NEARESTFACE:
{
if (samplePatch().size() && samplePatch() != "none")
{
FatalErrorIn
(
"mappedPatchBase::findSamples(const pointField&,"
" labelList&, labelList&, pointField&) const"
) << "No need to supply a patch name when in "
<< sampleModeNames_[mode] << " mode." << exit(FatalError);
//- Note: face-diagonal decomposition
const meshSearchMeshObject& meshSearchEngine =
meshSearchMeshObject::New(mesh);
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
forAll(samples, sampleI)
{
const point& sample = samples[sampleI];
label faceI = meshSearchEngine.findNearestFace(sample);
if (faceI == -1)
{
nearest[sampleI].second().first() = Foam::sqr(GREAT);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
else
{
const point& fc = mesh.faceCentres()[faceI];
nearest[sampleI].first() = pointIndexHit
(
true,
fc,
faceI
);
nearest[sampleI].second().first() = magSqr(fc-sample);
nearest[sampleI].second().second() = Pstream::myProcNo();
}
}
break;
}
case NEARESTPATCHFACEAMI:
{
// nothing to do here
return;
}
default:
{
FatalErrorIn("mappedPatchBase::findSamples(..)")
<< "problem." << abort(FatalError);
}
}
// Find nearest. Combine on master.
Pstream::listCombineGather(nearest, nearestEqOp());
Pstream::listCombineScatter(nearest);
Info<< "mappedPatchBase::findSamples on mesh " << sampleRegion()
<< " : " << endl;
forAll(nearest, sampleI)
{
label procI = nearest[sampleI].second().second();
label localI = nearest[sampleI].first().index();
Info<< " " << sampleI << " coord:"<< samples[sampleI]
<< " found on processor:" << procI
<< " in local cell/face/point:" << localI
<< " with location:" << nearest[sampleI].first().rawPoint()
<< endl;
}
}
// Convert back into proc+local index
sampleProcs.setSize(samples.size());
sampleIndices.setSize(samples.size());
sampleLocations.setSize(samples.size());
forAll(nearest, sampleI)
{
if (!nearest[sampleI].first().hit())
{
sampleProcs[sampleI] = -1;
sampleIndices[sampleI] = -1;
sampleLocations[sampleI] = vector::max;
}
else
{
sampleProcs[sampleI] = nearest[sampleI].second().second();
sampleIndices[sampleI] = nearest[sampleI].first().index();
sampleLocations[sampleI] = nearest[sampleI].first().hitPoint();
}
}
}
void Foam::mappedPatchBase::calcMapping() const
{
static bool hasWarned = false;
if (mapPtr_.valid())
{
FatalErrorIn("mappedPatchBase::calcMapping() const")
<< "Mapping already calculated" << exit(FatalError);
}
// Get points on face (since cannot use face-centres - might be off
// face-diagonal decomposed tets.
tmp<pointField> patchPoints(facePoints(patch_));
// Get offsetted points
const pointField offsettedPoints(samplePoints(patchPoints()));
andy
committed
// Do a sanity check - am I sampling my own patch?
// This only makes sense for a non-zero offset.
bool sampleMyself =
(
mode_ == NEARESTPATCHFACE
&& sampleRegion() == patch_.boundaryMesh().mesh().name()
&& samplePatch() == patch_.name()
);
// Check offset
vectorField d(offsettedPoints-patchPoints());
bool coincident = (gAverage(mag(d)) <= ROOTVSMALL);
if (sampleMyself && coincident)
{
WarningIn
(
"mappedPatchBase::mappedPatchBase\n"
"(\n"
" const polyPatch& pp,\n"
" const word& sampleRegion,\n"
" const sampleMode mode,\n"
" const word& samplePatch,\n"
" const vector& offset\n"
")\n"
) << "Invalid offset " << d << endl
<< "Offset is the vector added to the patch face centres to"
<< " find the patch face supplying the data." << endl
<< "Setting it to " << d
<< " on the same patch, on the same region"
<< " will find the faces themselves which does not make sense"
<< " for anything but testing." << endl
<< "patch_:" << patch_.name() << endl
<< "sampleRegion_:" << sampleRegion() << endl
<< "mode_:" << sampleModeNames_[mode_] << endl
<< "samplePatch_:" << samplePatch() << endl
<< "offsetMode_:" << offsetModeNames_[offsetMode_] << endl;
}
// Get global list of all samples and the processor and face they come from.
pointField samples;
labelList patchFaceProcs;
labelList patchFaces;
pointField patchFc;
collectSamples
(
patchPoints,
samples,
patchFaceProcs,
patchFaces,
patchFc
);
// Find processor and cell/face samples are in and actual location.
labelList sampleProcs;
labelList sampleIndices;
pointField sampleLocations;
findSamples(mode_, samples, sampleProcs, sampleIndices, sampleLocations);
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
// Check for samples that were not found. This will only happen for
// NEARESTCELL since finds cell containing a location
if (mode_ == NEARESTCELL)
{
label nNotFound = 0;
forAll(sampleProcs, sampleI)
{
if (sampleProcs[sampleI] == -1)
{
nNotFound++;
}
}
reduce(nNotFound, sumOp<label>());
if (nNotFound > 0)
{
if (!hasWarned)
{
WarningIn
(
"mappedPatchBase::mappedPatchBase\n"
"(\n"
" const polyPatch& pp,\n"
" const word& sampleRegion,\n"
" const sampleMode mode,\n"
" const word& samplePatch,\n"
" const vector& offset\n"
")\n"
) << "Did not find " << nNotFound
<< " out of " << sampleProcs.size() << " total samples."
<< " Sampling these on owner cell centre instead." << endl
<< "On patch " << patch_.name()
<< " on region " << sampleRegion()
<< " in mode " << sampleModeNames_[mode_] << endl
<< "with offset mode " << offsetModeNames_[offsetMode_]
<< ". Suppressing further warnings from " << type() << endl;
hasWarned = true;
}
// Collect the samples that cannot be found
DynamicList<label> subMap;
DynamicField<point> subSamples;
forAll(sampleProcs, sampleI)
{
if (sampleProcs[sampleI] == -1)
{
subMap.append(sampleI);
subSamples.append(samples[sampleI]);
// And re-search for pure nearest (should not fail)
labelList subSampleProcs;
labelList subSampleIndices;
pointField subSampleLocations;
findSamples
(
NEARESTONLYCELL,
subSamples,
subSampleProcs,
subSampleIndices,
subSampleLocations
);
// Insert
UIndirectList<label>(sampleProcs, subMap) = subSampleProcs;
UIndirectList<label>(sampleIndices, subMap) = subSampleIndices;
UIndirectList<point>(sampleLocations, subMap) = subSampleLocations;
// Now we have all the data we need:
// - where sample originates from (so destination when mapping):
// patchFaces, patchFaceProcs.
// - cell/face sample is in (so source when mapping)
// sampleIndices, sampleProcs.
//forAll(samples, i)
//{
// Info<< i << " need data in region "
// << patch_.boundaryMesh().mesh().name()
// << " for proc:" << patchFaceProcs[i]
// << " face:" << patchFaces[i]
// << " at:" << patchFc[i] << endl
// << "Found data in region " << sampleRegion()
// << " at proc:" << sampleProcs[i]
// << " face:" << sampleIndices[i]
// << " at:" << sampleLocations[i]
// << nl << endl;
//}
if (debug && Pstream::master())
{
OFstream str
(
patch_.boundaryMesh().mesh().time().path()
/ patch_.name()
+ "_mapped.obj"
);
Pout<< "Dumping mapping as lines from patch faceCentres to"
<< " sampled cell/faceCentres/points to file " << str.name()
<< endl;
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
label vertI = 0;
forAll(patchFc, i)
{
meshTools::writeOBJ(str, patchFc[i]);
vertI++;
meshTools::writeOBJ(str, sampleLocations[i]);
vertI++;
str << "l " << vertI-1 << ' ' << vertI << nl;
}
}
// Determine schedule.
mapPtr_.reset(new mapDistribute(sampleProcs, patchFaceProcs));
// Rework the schedule from indices into samples to cell data to send,
// face data to receive.
labelListList& subMap = mapPtr_().subMap();
labelListList& constructMap = mapPtr_().constructMap();
forAll(subMap, procI)
{
subMap[procI] = UIndirectList<label>
(
sampleIndices,
subMap[procI]
);
constructMap[procI] = UIndirectList<label>
(
patchFaces,
constructMap[procI]
);
//if (debug)
//{
// Pout<< "To proc:" << procI << " sending values of cells/faces:"
// << subMap[procI] << endl;
// Pout<< "From proc:" << procI
// << " receiving values of patch faces:"
// << constructMap[procI] << endl;
//}
}
// Redo constructSize
mapPtr_().constructSize() = patch_.size();
if (debug)
{
// Check that all elements get a value.
PackedBoolList used(patch_.size());
forAll(constructMap, procI)
{
const labelList& map = constructMap[procI];
forAll(map, i)
{
label faceI = map[i];
if (used[faceI] == 0)
{
used[faceI] = 1;
}
else
{
FatalErrorIn("mappedPatchBase::calcMapping() const")
<< "On patch " << patch_.name()
<< " patchface " << faceI
<< " is assigned to more than once."
<< abort(FatalError);
}
}
}
forAll(used, faceI)
{
if (used[faceI] == 0)
{
FatalErrorIn("mappedPatchBase::calcMapping() const")
<< "On patch " << patch_.name()
<< " patchface " << faceI
<< " is never assigned to."
<< abort(FatalError);
}
}
}
}
const Foam::autoPtr<Foam::searchableSurface>& Foam::mappedPatchBase::surfPtr()
const
{
const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
if (!surfPtr_.valid() && surfType != "none")
{
word surfName(surfDict_.lookupOrDefault("name", patch_.name()));
const polyMesh& mesh = patch_.boundaryMesh().mesh();
surfPtr_ =
searchableSurface::New
(
surfType,
IOobject
(
surfName,
mesh.time().constant(),
"triSurface",
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
),
surfDict_
);
}
return surfPtr_;
}
void Foam::mappedPatchBase::calcAMI() const
{
if (AMIPtr_.valid())
{
FatalErrorIn("mappedPatchBase::calcAMI() const")
<< "AMI already calculated" << exit(FatalError);
}
AMIPtr_.clear();
andy
committed
const polyPatch& nbr = samplePolyPatch();
andy
committed
// Transform neighbour patch to local system
pointField nbrPoints(samplePoints(nbr.localPoints()));
primitivePatch nbrPatch0
(
SubList<face>
nbr.localFaces(),
nbr.size()
),
nbrPoints
);
if (debug)
{
OFstream os(patch_.name() + "_neighbourPatch-org.obj");
meshTools::writeOBJ(os, samplePolyPatch().localFaces(), nbrPoints);
OFstream osN(patch_.name() + "_neighbourPatch-trans.obj");
meshTools::writeOBJ(osN, nbrPatch0, nbrPoints);
OFstream osO(patch_.name() + "_ownerPatch.obj");
meshTools::writeOBJ(osO, patch_.localFaces(), patch_.localPoints());
}
andy
committed
// Construct/apply AMI interpolation to determine addressing and weights
AMIPtr_.reset
(
new AMIPatchToPatchInterpolation
(
patch_,
faceAreaIntersect::tmMesh,
)
);
}
mattijs
committed
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
// Hack to read old (List-based) format. See Field.C. The difference
// is only that in case of falling back to old format it expects a non-uniform
// list instead of a single vector.
Foam::tmp<Foam::pointField> Foam::mappedPatchBase::readListOrField
(
const word& keyword,
const dictionary& dict,
const label size
)
{
tmp<pointField> tfld(new pointField());
pointField& fld = tfld();
if (size)
{
ITstream& is = dict.lookup(keyword);
// Read first token
token firstToken(is);
if (firstToken.isWord())
{
if (firstToken.wordToken() == "uniform")
{
fld.setSize(size);
fld = pTraits<vector>(is);
}
else if (firstToken.wordToken() == "nonuniform")
{
is >> static_cast<List<vector>&>(fld);
if (fld.size() != size)
{
FatalIOErrorIn
(
"mappedPatchBase::readListOrField"
"(const word& keyword, const dictionary&, const label)",
dict
) << "size " << fld.size()
<< " is not equal to the given value of " << size
<< exit(FatalIOError);
}
}
else
{
FatalIOErrorIn
(
"mappedPatchBase::readListOrField"
"(const word& keyword, const dictionary&, const label)",
dict
) << "expected keyword 'uniform' or 'nonuniform', found "
<< firstToken.wordToken()
<< exit(FatalIOError);
}
}
else
{
if (is.version() == 2.0)
{
IOWarningIn
(
"mappedPatchBase::readListOrField"
"(const word& keyword, const dictionary&, const label)",
dict
) << "expected keyword 'uniform' or 'nonuniform', "
"assuming List format for backwards compatibility."
"Foam version 2.0." << endl;
is.putBack(firstToken);
is >> static_cast<List<vector>&>(fld);
}
}
}
return tfld;
}
// * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
Foam::mappedPatchBase::mappedPatchBase
(
const polyPatch& pp
)
:
patch_(pp),
sampleRegion_(patch_.boundaryMesh().mesh().name()),
mode_(NEARESTPATCHFACE),
samplePatch_(""),
coupleGroup_(),
offsetMode_(UNIFORM),
offset_(vector::zero),
offsets_(pp.size(), offset_),
distance_(0),
sameRegion_(sampleRegion_ == patch_.boundaryMesh().mesh().name()),
mapPtr_(NULL),
AMIPtr_(NULL),
AMIReverse_(false),
surfDict_(fileName("surface"))
{}
Foam::mappedPatchBase::mappedPatchBase
(
const polyPatch& pp,
const word& sampleRegion,
const sampleMode mode,
const word& samplePatch,
const vectorField& offsets
)
:
patch_(pp),
sampleRegion_(sampleRegion),
mode_(mode),
samplePatch_(samplePatch),