adjointRotatingWallVelocityFvPatchVectorField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020 PCOpt/NTUA
9 Copyright (C) 2020 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34Foam::adjointRotatingWallVelocityFvPatchVectorField::
35adjointRotatingWallVelocityFvPatchVectorField
36(
37 const fvPatch& p,
39)
40:
42 origin_(),
43 axis_(Zero),
44 omega_(nullptr)
45{}
46
47
48Foam::adjointRotatingWallVelocityFvPatchVectorField::
49adjointRotatingWallVelocityFvPatchVectorField
50(
52 const fvPatch& p,
54 const fvPatchFieldMapper& mapper
55)
56:
58 origin_(ptf.origin_),
59 axis_(ptf.axis_),
60 omega_(ptf.omega_.clone())
61{}
62
63
64Foam::adjointRotatingWallVelocityFvPatchVectorField::
65adjointRotatingWallVelocityFvPatchVectorField
66(
67 const fvPatch& p,
69 const dictionary& dict
70)
71:
73 origin_(dict.get<vector>("origin")),
74 axis_(dict.get<vector>("axis")),
75 omega_(Function1<scalar>::New("omega", dict, &db()))
76{}
77
78
79Foam::adjointRotatingWallVelocityFvPatchVectorField::
80adjointRotatingWallVelocityFvPatchVectorField
81(
84)
85:
87 origin_(pivpvf.origin_),
88 axis_(pivpvf.axis_),
89 omega_(pivpvf.omega_.clone())
90{}
91
92
93// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94
97{
98 const scalar t(this->db().time().timeOutputValue());
99 const scalar om(omega_->value(t));
100 const vector omega(om*axis_/mag(axis_));
101 tensor mult
102 (
103 scalar(0), -omega.z(), omega.y(),
104 omega.z(), scalar(0), -omega.x(),
105 -omega.y(), omega.x(), scalar(0)
106 );
107
108 return tmp<tensorField>::New(patch().size(), mult);
109}
110
111
113(
114 Ostream& os
115) const
116{
118 os.writeEntry("origin", origin_);
119 os.writeEntry("axis", axis_);
120 omega_->writeData(os);
121}
122
123
124// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125
126namespace Foam
127{
129 (
132 );
133}
134
135// ************************************************************************* //
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...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
The same as adjointWallVelocity but additionally computes the sensitivity contribution emerging from ...
virtual tmp< tensorField > dxdbMult() const
Compute contribution to SDs.
Adjoint wall velocity boundary condition. If nutUSpaldingWallFunction is employed in the flow solutio...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dictionary dict
optimisationManager & om
Definition: createFields.H:6