Go to the documentation of this file.
44 fixedValueFvPatchVectorField(
p, iF),
60 fixedValueFvPatchVectorField(ptf,
p, iF, mapper),
75 fixedValueFvPatchVectorField(
p, iF),
77 kappa_(
dict.getOrDefault<scalar>(
"kappa", 0.41)),
78 E_(
dict.getOrDefault<scalar>(
"E", 9.8))
94 fixedValueFvPatchVectorField(pivpvf, iF),
96 kappa_(pivpvf.kappa_),
113 SAwallFunctionPatchField;
117 isA<SAwallFunctionPatchField>(boundaryContrPtr_->turbulentDiffusivity())
118 &&
patch().size() != 0
147 const scalarField auxA((kappa_/E_)*(
exp(kUu)-1 - kUu - 0.5*kUu*kUu));
157 label cellI =
patch().faceCells()[faceI];
159 2*auxB[faceI]*vtau[faceI]*((rt[faceI] + Uap_t[faceI]))
160 *(Uc[faceI]/
mag(Uc[faceI]))*magSf[faceI];
173 const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
176 tmp<vectorField> tnf =
patch().nf();
191 tmp<vectorField> tsource = boundaryContrPtr_->normalVelocitySource();
198 SAwallFunctionPatchField;
201 if (isA<SAwallFunctionPatchField>(nutb))
203 Uap_t1 = (Uac & tf1)*tf1;
209 Uap_t1 = - (tsource() & tf1)*tf1;
217 fixedValueFvPatchVectorField::updateCoeffs();
224 writeEntry(
"value",
os);
238 adjointWallVelocityFvPatchVectorField
adjointWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
virtual tmp< Field< Type > > snGrad() const
A class for handling words, derived from Foam::string.
virtual void write(Ostream &) const
A class for managing temporary objects.
static constexpr const zero Zero
Base class for solution control classes.
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
dimensionedScalar exp(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
virtual void manipulateMatrix(fvMatrix< vector > &matrix)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Adjoint wall velocity boundary condition. If nutUSpaldingWallFunction is employed in the flow solutio...
virtual tmp< Field< Type > > patchInternalField() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
adjointBoundaryCondition< vector > adjointVectorBoundaryCondition
fvPatchField< vector > fvPatchVectorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Ostream & writeEntry(const keyType &key, const T &value)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...