37Foam::turbulentDigitalFilterInletFvPatchField<Type>::calcPatchNormal()
const
45 const scalar error =
max(
magSqr(patchNormal + nf));
50 <<
"Patch " << this->
patch().name() <<
" is not planar"
54 return patchNormal.normalise();
59void Foam::turbulentDigitalFilterInletFvPatchField<Type>::initialisePatch()
63 AMIPtr_->calculate(this->
patch().
patch(), L_.patch());
65 patchNormal_ = calcPatchNormal();
70void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapL
75 Field<Type> sourceFld;
77 if (Pstream::master())
79 sourceFld = L_.convolve();
86 AMIPtr_->interpolateToSource
89 multiplyWeightedOp<Type, plusEqOp<Type>>(cop),
103void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapR
108 const scalar t = this->db().time().timeOutputValue();
120void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapR
125 const scalar t = this->db().time().timeOutputValue();
129 for (symmTensor& r :
R)
136 r.yz() = (r.yz() - r.xy()*r.xz())/r.yy();
147 u.z() = u.x()*r.xz() + u.y()*r.yz() + u.z()*r.zz();
148 u.y() = u.x()*r.xy() + u.y()*r.yy();
149 u.x() = u.x()*r.xx();
155void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapMean
160 const scalar t = this->db().time().timeOutputValue();
162 fld += meanPtr_->value(t);
167void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapMean
172 const scalar t = this->db().time().timeOutputValue();
173 tmp<vectorField> tmean = meanPtr_->value(t);
181 gSum((bulk & patchNormal_)*this->
patch().magSf())
224 AMIPtr_(ptf.AMIPtr_.clone()),
225 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
226 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
227 curTimeIndex_(ptf.curTimeIndex_),
228 patchNormal_(ptf.patchNormal_),
256 this->patch().patch(),
265 this->patch().patch(),
277 if (!L_.fsm() && this->db().time().isAdjustTimeStep())
280 <<
"Varying time-step computations are not "
281 <<
"supported by the digital filter method."
285 const scalar t = this->db().time().timeOutputValue();
300 AMIPtr_(ptf.AMIPtr_.clone()),
301 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
302 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
303 curTimeIndex_(ptf.curTimeIndex_),
304 patchNormal_(ptf.patchNormal_),
318 AMIPtr_(ptf.AMIPtr_.clone()),
319 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
320 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
321 curTimeIndex_(ptf.curTimeIndex_),
322 patchNormal_(ptf.patchNormal_),
339 meanPtr_->autoMap(m);
351 const fvPatchField<Type>& ptf,
362 meanPtr_->rmap(dfmptf.meanPtr_(), addr);
366 Rptr_->rmap(dfmptf.Rptr_(), addr);
379 if (curTimeIndex_ == -1)
384 if (curTimeIndex_ != this->db().time().
timeIndex())
395 curTimeIndex_ = this->db().time().timeIndex();
412 meanPtr_->writeData(
os);
416 Rptr_->writeData(
os);
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Face area weighted Arbitrary Mesh Interface (AMI) method.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
static int debug
Flag to activate debug statements.
static void checkStresses(const symmTensorField &R)
Check if input Reynolds stresses are valid.
Digital-filter based boundary condition for vector- and scalar-based quantities (e....
virtual void write(Ostream &) const
Write boundary condition settings.
virtual void rmap(const fvPatchField< Type > &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
turbulentDigitalFilterInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
static constexpr const zero Zero
Global zero (0).
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
#define forAll(list, i)
Loop across all elements in list.