64#define SNAP_END_ENCODE(pos, val) (((val) << (4 + 2 * pos)))
65#define SNAP_END_VALUE(pos, val) (((val) >> (4 + 2 * pos)) & 0x3)
93 if (doSnap && cutIndex && cutIndex != 0xF)
98 p0 -= val;
if (cutIndex & 1)
p0 *= -1;
99 p1 -= val;
if (cutIndex & 2) p1 *= -1;
100 p2 -= val;
if (cutIndex & 4) p2 *= -1;
101 p3 -= val;
if (cutIndex & 8) p3 *= -1;
106 #undef ADD_SNAP_INDEX
107 #define ADD_SNAP_INDEX(pos, d1, d2, idx1, idx2) \
108 switch (cutIndex & (idx1 | idx2)) \
112 cutIndex |= SNAP_END_ENCODE \
115 ((d1 * 100 < d2) ? 1 : (d2 * 100 < d1) ? 2 : 0) \
126 #undef ADD_SNAP_INDEX
137 DynamicList<label>& verts,
144 if (a !=
b &&
b != c && c != a)
181Foam::isoSurfaceTopo::tetCutAddressing::tetCutAddressing
183 const label nCutCells,
185 const bool useDebugCuts
188 vertsToPointLookup_(12*nCutCells),
191 pointToFace_(10*nCutCells),
192 pointFromDiag_(10*nCutCells),
194 pointToVerts_(10*nCutCells),
195 cutPoints_(12*nCutCells),
198 debugCutTetsOn_(useDebugCuts)
209 snapVertsLookup_.reserve(2*nCutCells);
213 debugCutTets_.reserve(6*nCutCells);
220void Foam::isoSurfaceTopo::tetCutAddressing::clearDebug()
222 debugCutTets_.clearStorage();
226void Foam::isoSurfaceTopo::tetCutAddressing::clearDiagonal()
228 pointToFace_.clearStorage();
229 pointFromDiag_.clearStorage();
233void Foam::isoSurfaceTopo::tetCutAddressing::clearHashes()
235 vertsToPointLookup_.clear();
236 snapVertsLookup_.clear();
240Foam::label Foam::isoSurfaceTopo::tetCutAddressing::generatePoint
251 label pointi = vertsToPointLookup_.lookup(
vertices, -1);
254 bool addNewPoint(
true);
256 const label snapPointi =
259 : (snapEnd == 2) ?
vertices.second()
263 if (snapPointi == -1)
266 pointi = pointToVerts_.size();
273 edgeIsDiagonal =
false;
275 pointi = snapVertsLookup_.lookup(snapPointi, -1);
276 addNewPoint = (pointi == -1);
279 pointi = pointToVerts_.size();
280 snapVertsLookup_.insert(snapPointi, pointi);
281 pointToVerts_.append(edge(snapPointi, snapPointi));
287 pointToFace_.append(facei);
288 pointFromDiag_.append(edgeIsDiagonal);
291 vertsToPointLookup_.insert(
vertices, pointi);
298bool Foam::isoSurfaceTopo::tetCutAddressing::generatePoints
301 const int tetCutIndex,
302 const tetCell& tetLabels,
305 const FixedList<bool, 6>& edgeIsDiagonal
309 const label nCutPointsOld(cutPoints_.size());
312 switch (tetCutIndex & 0x0F)
319 case 0x0E: flip =
true; [[fallthrough]];
358 case 0x0D: flip =
true; [[fallthrough]];
397 case 0x0C: flip =
true; [[fallthrough]];
447 case 0x0B: flip =
true; [[fallthrough]];
486 case 0x0A: flip =
true; [[fallthrough]];
536 case 0x09: flip =
true; [[fallthrough]];
586 case 0x07: flip =
true; [[fallthrough]];
625 const bool added(nCutPointsOld != cutPoints_.size());
627 if (added && debugCutTetsOn_)
629 debugCutTets_.append(tetLabels.shape());
639void Foam::isoSurfaceTopo::generateTriPoints
644 tetCutAddressing& tetCutAddr
649 const cell& cFaces =
mesh_.cells()[celli];
650 const bool doSnap = this->
snap();
657 const label facei = cFaces[0];
658 const face& f0 = faces[facei];
662 const face& f1 = faces[cFaces[1]];
667 if (!f0.found(apexi))
673 const label
p0 = f0[0];
677 if (faceOwner[facei] == celli)
682 const tetCell tetLabels(
p0, p1, p2, apexi);
683 const int tetCutIndex
696 tetCutAddr.generatePoints
701 FixedList<bool, 6>(
false)
706 for (
const label facei : cFaces)
710 !
mesh_.isInternalFace(facei)
717 const face&
f = faces[facei];
719 label fp0 = tetBasePtIs[facei];
727 const label
p0 =
f[fp0];
728 label fp =
f.fcIndex(fp0);
729 for (label i = 2; i <
f.size(); ++i)
735 FixedList<bool, 6> edgeIsDiagonal(
false);
736 if (faceOwner[facei] == celli)
739 if (i != 2) edgeIsDiagonal[1] =
true;
740 if (i !=
f.size()-1) edgeIsDiagonal[0] =
true;
744 if (i != 2) edgeIsDiagonal[0] =
true;
745 if (i !=
f.size()-1) edgeIsDiagonal[1] =
true;
748 const tetCell tetLabels(
p0, p1, p2,
mesh_.nPoints()+celli);
749 const int tetCutIndex
762 tetCutAddr.generatePoints
777void Foam::isoSurfaceTopo::triangulateOutside
779 const bool filterDiag,
786 DynamicList<face>& compactFaces,
787 DynamicList<label>& compactCellIDs
804 face&
f = compactFaces.emplace_back(loop.size());
809 const label pointi =
mp[loop[i]];
810 if (filterDiag && pointFromDiag[pointi])
812 const label prevPointi =
mp[loop[loop.fcIndex(i)]];
815 pointFromDiag[prevPointi]
816 && (pointToFace[pointi] != pointToFace[prevPointi])
841 const label pointi =
mp[loop[i]];
845 compactCellIDs.append(cellID);
851void Foam::isoSurfaceTopo::removeInsidePoints
854 const bool filterDiag,
862 DynamicList<label>& compactCellIDs
867 compactCellIDs.clear();
868 compactCellIDs.reserve(
s.size()/4);
870 DynamicList<face> compactFaces(
s.size()/4);
872 for (label celli = 0; celli < start.size()-1; ++celli)
876 const label nTris = start[celli+1]-start[celli];
882 SubList<face>(
s, nTris, start[celli]),
900 s.swapFaces(compactFaces);
926 Pout<<
"isoSurfaceTopo:" <<
nl
929 <<
" isoValue : " << iso <<
nl
932 <<
" mesh span : " <<
mesh.bounds().mag() <<
nl
933 <<
" ignoreCells : " << ignoreCells.
count()
940 label nBlockedCells = 0;
943 nBlockedCells +=
blockCells(cellCutType_, ignoreCells);
961 Pout<<
"isoSurfaceTopo : candidate cells cut "
963 <<
" blocked " << nBlockedCells
975 "isoSurfaceTopo.cutType",
976 fvmesh.time().timeName(),
986 auto& debugFld = debugField.primitiveFieldRef();
988 forAll(cellCutType_, celli)
990 debugFld[celli] = cellCutType_[celli];
993 Info<<
"Writing cut types: " << debugField.objectRelPath() <<
nl;
1006 forAll(cellCutType_, celli)
1008 if ((cellCutType_[celli] & realCut) != 0)
1010 debugCutCells[nout] = celli;
1012 if (nout >= nCutCells)
break;
1017 vtk::vtuCells vtuCells;
1018 vtuCells.reset(
mesh_, debugCutCells);
1020 vtk::internalMeshWriter
writer
1026 mesh_.time().globalPath()
1027 / (
"isoSurfaceTopo." + timeDesc +
"-cutCells")
1042 Info<<
"isoSurfaceTopo : (debug) wrote "
1050 tetCutAddressing tetCutAddr
1058 for (label celli = 0; celli <
mesh_.nCells(); ++celli)
1060 startTri[celli] = tetCutAddr.nFaces();
1074 startTri.back() = tetCutAddr.nFaces();
1077 tetBasePtIs.clear();
1078 tetCutAddr.clearHashes();
1082 faceList allTriFaces(startTri.back());
1084 auto& verts = tetCutAddr.cutPoints();
1087 for (face&
f : allTriFaces)
1090 f[0] = verts[verti++];
1091 f[1] = verts[verti++];
1092 f[2] = verts[verti++];
1094 verts.clearStorage();
1100 for (label celli = 0; celli < startTri.size()-1; ++celli)
1106 (startTri[celli+1] - startTri[celli]),
1112 pointToVerts_.transfer(tetCutAddr.pointToVerts());
1118 mesh_.cellCentres(),
1127 std::move(allTriPoints),
1128 std::move(allTriFaces),
1134 Pout<<
"isoSurfaceTopo : generated "
1140 if ((debug & 8) && (params.
filter() != filterType::NONE))
1142 const Mesh&
s = *
this;
1144 vtk::surfaceWriter
writer
1150 mesh_.time().globalPath()
1151 / (
"isoSurfaceTopo." + timeDesc +
"-triangles")
1172 const edge& verts = pointToVerts_[i];
1173 if (verts.first() == verts.second())
1178 if (tetCutAddr.pointFromDiag().test(i))
1181 pointStatus[i] = -1;
1185 writer.write(
"point-status", pointStatus);
1188 Info<<
"isoSurfaceTopo : (debug) wrote "
1203 if (params.
filter() == filterType::NONE)
1211 s.compactPoints(pointMap);
1212 pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1223 DynamicList<label> compactCellIDs;
1230 params.
filter() == filterType::DIAGCELL
1231 || params.
filter() == filterType::NONMANIFOLD
1233 tetCutAddr.pointFromDiag(),
1234 tetCutAddr.pointToFace(),
1240 s.compactPoints(pointMap);
1242 pointToVerts_ = UIndirectList<edge>(pointToVerts_, pointMap)();
1247 Pout<<
"isoSurfaceTopo :"
1248 " after removing cell centre and face-diag triangles: "
1256 tetCutAddr.clearDiagonal();
1262 bitSet isProtectedPoint;
1265 (params.
filter() == filterType::NONMANIFOLD)
1266 || tetCutAddr.debugCutTetsOn()
1272 isProtectedPoint.resize(
mesh_.nPoints());
1276 label facei =
mesh_.nInternalFaces();
1277 facei <
mesh_.nFaces();
1281 isProtectedPoint.set(
mesh_.faces()[facei]);
1290 for (label facei = 0; facei <
mesh.nInternalFaces(); ++facei)
1300 isProtectedPoint.set(
mesh_.faces()[facei]);
1307 cellCutType_.
clear();
1311 if (tetCutAddr.debugCutTetsOn())
1315 const auto& debugCuts = tetCutAddr.debugCutTets();
1318 vtk::vtuCells vtuCells;
1319 vtuCells.resetShapes(debugCuts);
1322 vtuCells.setNumPoints(
mesh_.nPoints());
1325 vtk::internalMeshWriter
writer
1331 mesh_.time().globalPath()
1332 / (
"isoSurfaceTopo." + timeDesc +
"-cutTets")
1344 Field<scalar> cutTetQuality(debugCuts.size());
1345 forAll(cutTetQuality, teti)
1355 writer.writeCellData(
"tetQuality", cutTetQuality);
1365 for (
const edge& verts : pointToVerts_)
1367 if (verts.first() == verts.second())
1370 pointStatus[verts.first()] = 1;
1374 writer.writePointData(
"point-status", pointStatus);
1377 Info<<
"isoSurfaceTopo : (debug) wrote "
1388 const Mesh&
s = *
this;
1393 openEdgeIds.reserve(
s.size());
1398 if (eFaces.size() == 1)
1402 const edge&
e = surfEdges[edgei];
1403 const edge& verts0 = pointToVerts_[
mp[
e.first()]];
1404 const edge& verts1 = pointToVerts_[
mp[
e.second()]];
1408 isProtectedPoint.test(verts0.first())
1409 && isProtectedPoint.test(verts0.second())
1410 && isProtectedPoint.test(verts1.first())
1411 && isProtectedPoint.test(verts1.second())
1419 openEdgeIds.insert(edgei);
1424 const label nOpenEdges
1434 openEdgeIds.sortedToc()
1443 mesh_.time().globalPath()
1444 / (
"isoSurfaceTopo." + timeDesc +
"-openEdges")
1454 Info<<
"isoSurfaceTopo : (debug) wrote "
1455 << nOpenEdges <<
" open edges: "
1460 Info<<
"isoSurfaceTopo : no open edges" <<
nl;
1467 const Mesh&
s = *
this;
1469 vtk::surfaceWriter
writer
1475 mesh_.time().globalPath()
1476 / (
"isoSurfaceTopo." + timeDesc +
"-surface")
1497 const edge& verts = pointToVerts_[i];
1498 if (verts.first() == verts.second())
1505 writer.write(
"point-status", pointStatus);
1508 Info<<
"isoSurfaceTopo : (debug) wrote "
1514 tetCutAddr.clearDebug();
1517 if (params.
filter() == filterType::NONMANIFOLD)
1531 bitSet faceSelection;
1538 UIndirectList<face>(surf, faceAddr),
1542 faceSelection.clear();
1543 faceSelection.resize(erosion.size());
1546 const edgeList& surfEdges = erosion.edges();
1549 label nEdgeRemove = 0;
1554 if (eFaces.size() == 1)
1558 const edge&
e = surfEdges[edgei];
1559 const edge& verts0 = pointToVerts_[
mp[
e.first()]];
1560 const edge& verts1 = pointToVerts_[
mp[
e.second()]];
1564 isProtectedPoint.test(verts0.first())
1565 && isProtectedPoint.test(verts0.second())
1566 && isProtectedPoint.test(verts1.first())
1567 && isProtectedPoint.test(verts1.second())
1575 faceSelection.set(eFaces[0]);
1583 Pout<<
"isoSurfaceTopo :"
1584 <<
" removing " << faceSelection.count()
1585 <<
" / " << faceSelection.size()
1586 <<
" faces on " << nEdgeRemove <<
" open edges" <<
endl;
1602 if (surf.size() != faceAddr.size())
1604 faceSelection.clear();
1605 faceSelection.resize(surf.size());
1606 faceSelection.set(faceAddr);
1622 Mesh::subsetMesh(include, pointMap,
faceMap)
1624 Mesh::transfer(filtered);
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
void clearStorage()
Clear the list and delete storage.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
label size() const noexcept
The number of elements in table.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
fileName objectRelPath() const
The object path relative to the case.
label size() const noexcept
The number of elements in the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void clear()
Clear the list, i.e. set size to zero.
MeshedSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
void transfer(pointField &pointLst, List< face > &faceLst)
void resize(const label numElem, const unsigned int val=0u)
Reset addressable list size, does not shrink the allocated size.
const T & first() const noexcept
Access the first element.
const T & second() const noexcept
Access the second element.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
const labelListList & edgeLoops() const
const labelListList & edgeFaces() const
label timeIndex() const noexcept
Return the current time index.
A List with indirect addressing. Like IndirectList but does not store addressing.
T & first()
Access first element of the list, position [0].
bool test(const label i) const
Test bool value at specified position, always false for out-of-range access.
T & back()
Access last element of the list, position [size()-1].
void size(const label n)
Older name for setAddressableSize.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
unsigned int count(const bool on=true) const
Count number of bits set.
void set(const bitSet &bitset)
Set specified bits from another bitset.
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Face selection method for createBaffles.
A face is a list of labels corresponding to mesh vertices.
A class for handling file names.
Low-level components common to various iso-surface algorithms.
bitSet ignoreBoundaryFaces_
Optional boundary faces to ignore.
isoSurfaceBase(const isoSurfaceBase &)=delete
No copy construct.
@ BLOCKED
Blocked (never cut).
@ ANYCUT
Any cut type (bitmask).
@ TETCUT
Cell cut is a tet.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
const scalarField & pVals_
Point values.
label calcCellCuts(List< cutType > &cuts) const
Populate a list of candidate cell cuts using getCellCutType().
labelList meshCells_
For every face, the original cell in mesh.
const scalar iso_
Iso value.
void ignoreCyclics()
Set ignoreBoundaryFaces to ignore cyclics (cyclicACMI).
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
const polyMesh & mesh_
Reference to mesh.
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
const scalarField & cVals_
Cell values.
Preferences for controlling iso-surface algorithms.
filterType filter() const noexcept
Get current filter type.
static const Enum< filterType > filterNames
Names for the filtering types.
const boundBox & getClipBounds() const noexcept
Get optional clipping bounding box.
bool snap() const noexcept
Get point snapping flag.
Marching tet iso surface algorithm with optional filtering to keep only points originating from mesh ...
void inplaceSubsetMesh(const bitSet &include)
Subset the surface using the selected faces.
isoSurfaceTopo(const polyMesh &mesh, const scalarField &cellValues, const scalarField &pointValues, const scalar iso, const isoSurfaceParams ¶ms=isoSurfaceParams(), const bitSet &ignoreCells=bitSet())
Construct from cell and point values.
tmp< Field< Type > > interpolateTemplate(const Field< Type > &cellData, const Field< Type > &pointData) const
Interpolates cellData and pointData fields.
const Time & time() const noexcept
Return time registry.
static labelList adjustTetBasePtIs(const polyMesh &mesh, const bool report=false)
Return an adjusted list of tet base points.
Mesh consisting of general polyhedral cells.
virtual const pointField & points() const
Return raw points.
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
scalar quality() const
Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron h...
@ OUTSIDE
A location outside the volume.
Write an OpenFOAM volume (internal) geometry and internal fields as a vtu file or a legacy vtk file.
virtual bool beginPointData(label nFields=0)
Begin PointData for specified number of fields.
Write edge/points (optionally with fields) as a vtp file or a legacy vtk file.
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
A deep-copy description of an OpenFOAM volume mesh in data structures suitable for VTK UnstructuredGr...
void resetShapes(const UList< cellShape > &shapes)
Reset sizing using primitive shapes only (ADVANCED USAGE).
void addPointCellLabels(const labelUList &cellIds)
Define which additional cell-centres are to be used (ADVANCED!).
void reset(const polyMesh &mesh)
Create the geometry using the previously requested output and decomposition types.
void setNumPoints(label n) noexcept
Alter number of mesh points (ADVANCED USAGE).
label nFieldPoints() const noexcept
Number of field points = nPoints + nAddPoints.
A class for handling words, derived from Foam::string.
static word printf(const char *fmt, const PrimitiveType &val)
Use a printf-style formatter for a primitive.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & p0
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Convenience macros for instantiating iso-surface interpolate methods.
#define defineIsoSurfaceInterpolateMethods(ThisClass)
#define ADD_SNAP_INDEX(pos, d1, d2, idx1, idx2)
#define SNAP_END_VALUE(pos, val)
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
List< edge > edgeList
List of edge.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
void inplaceSubset(const BoolListType &select, ListType &input, const bool invert=false)
Inplace extract elements of the input list when select is true.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
List< surfZone > surfZoneList
List of surfZone.
static void appendTriLabels(DynamicList< label > &verts, const label a, const label b, const label c, const bool flip)
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
List< face > faceList
List of faces.
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
pointField vertices(const blockVertexList &bvl)
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
UList< bool > boolUList
A UList of bools.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const point & getMeshPointRef(const polyMesh &mesh, const label pointi)
UList< label > labelUList
A UList of labels.
static int getTetCutIndex(scalar p0, scalar p1, scalar p2, scalar p3, const scalar val, const bool doSnap) noexcept
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.