38void Foam::pointAttractionDisplacementPointPatchVectorField::calcProjection
40 vectorField& displacement
52 const scalar projectLen =
mesh.bounds().mag();
57 const pointZone* zonePtr =
nullptr;
59 if (frozenPointsZone_.size() > 0)
63 zonePtr = &
pZones[frozenPointsZone_];
65 Pout<<
"pointAttractionDisplacementPointPatchVectorField : Fixing all "
66 << zonePtr->size() <<
" points in pointZone " << zonePtr->
name()
80 start[i] =
points0[meshPoints[i]] + displacement[i];
85 label nNotProjected = 0;
88 const label meshPointi = meshPoints[i];
89 const point& pt =
mesh.points()[meshPointi];
91 if (zonePtr && (zonePtr->whichPoint(meshPointi) >= 0))
94 displacement[i] =
points0[meshPointi] - pt;
101 displacement[i] = nearest.point() -
points0[meshPointi];
109 Pout<<
" point:" << meshPointi
111 <<
" did not find any surface within " << projectLen
118 reduce(nNotProjected, sumOp<label>());
120 if (nNotProjected > 0)
122 Info<<
"pointAttractionDisplacement :"
124 <<
" did not project " << nNotProjected
125 <<
" out of " <<
returnReduce(meshPoints.size(), sumOp<label>())
126 <<
" points." <<
endl;
156 frozenPointsZone_(
dict.getOrDefault(
"frozenPointsZone",
word::null))
177 velocity_(ppf.velocity_),
178 featFileName_(ppf.featFileName_),
179 frozenPointsZone_(ppf.frozenPointsZone_)
190 velocity_(ppf.velocity_),
191 featFileName_(ppf.featFileName_),
192 frozenPointsZone_(ppf.frozenPointsZone_)
204 velocity_(ppf.velocity_),
205 featFileName_(ppf.featFileName_),
206 frozenPointsZone_(ppf.frozenPointsZone_)
217 const Time& tm = this->
patch().boundaryMesh().mesh().time();
246 return pointTreePtr_();
257 const vectorField currentDisplacement(this->patchInternalField());
261 calcProjection(displacement);
265 vectorField offset(displacement-currentDisplacement);
269 const Time& tm = this->patch().boundaryMesh().mesh().time();
270 const scalar deltaT = tm.deltaTValue();
271 const vector clipVelocity = velocity_*deltaT;
277 const scalar magD(
mag(d));
278 if (magD > ROOTVSMALL)
281 d *=
min(magD,
mag(clipVelocity));
288 setInInternalField(iF, (currentDisplacement+offset)());
300 os.writeEntry(
"file", featFileName_);
301 os.writeEntryIfDifferent<
word>
307 os.writeEntry(
"velocity", velocity_);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
scalar deltaTValue() const noexcept
Return time step value.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Area discretisation.
static void read(const objectRegistry &obr, const dictionary &dict)
Read (& store) geometry. Exposed so point attraction can reuse it.
A class for handling file names.
Non-pointer based hierarchical recursive searching.
A class for handling keywords in dictionaries.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const word & name() const noexcept
The patch name.
Displacement by attraction to nearest point. Use in a displacementMotionSolver as a bc on the pointDi...
virtual void write(Ostream &) const
Write.
pointAttractionDisplacementPointPatchVectorField(const pointPatch &, const DimensionedField< vector, pointMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const indexedOctree< treeDataPoint > & pointTree() const
const pointMesh & mesh() const noexcept
Return the mesh reference.
const pointPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
Foam::pointPatchFieldMapper.
void patchInternalField(const UList< Type1 > &internalData, const labelUList &addressing, UList< Type1 > &pfld) const
void setInInternalField(Field< Type1 > &iF, const Field< Type1 > &pF, const labelUList &meshPoints) const
const Field< vector > & primitiveField() const noexcept
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Basic pointPatch represents a set of points from the mesh.
virtual const labelList & meshPoints() const =0
Return mesh points.
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Standard boundBox with extra functionality for use in octree.
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
static const word null
An empty word.
IOporosityModelList pZones(mesh)
OBJstream os(runTime.globalPath()/outputName)
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere).
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
pointPatchField< vector > pointPatchVectorField
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete pointPatchField type and add to run-time tables Example, (pointPatchScalarField,...
Tree tree(triangles.begin(), triangles.end())
#define forAll(list, i)
Loop across all elements in list.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))