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 -------------------------------------------------------------------------------
11 License
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 
36 (
37  const fvPatch& p,
39 )
40 :
42  origin_(),
43  axis_(Zero),
44  omega_(nullptr)
45 {}
46 
47 
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  adjointWallVelocityFvPatchVectorField(ptf, p, iF, mapper),
58  origin_(ptf.origin_),
59  axis_(ptf.axis_),
60  omega_(ptf.omega_.clone())
61 {}
62 
63 
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 
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 
126 namespace Foam
127 {
129  (
131  adjointRotatingWallVelocityFvPatchVectorField
132  );
133 }
134 
135 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Adjoint wall velocity boundary condition. If nutUSpaldingWallFunction is employed in the flow solutio...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
fvPatchField< vector > fvPatchVectorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
optimisationManager & om
Definition: createFields.H:6
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
The same as adjointWallVelocity but additionally computes the sensitivity contribution emerging from ...
Macros for easy insertion into run-time selection tables.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
adjointRotatingWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< tensorField > dxdbMult() const
Compute contribution to SDs.
Tensor of scalars, i.e. Tensor<scalar>.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127