Loading...
Searching...
No Matches
adjointLaminar.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) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2023 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "adjointLaminar.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
39namespace adjointRASModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44defineTypeNameAndDebug(adjointLaminar, 0);
45addToRunTimeSelectionTable(adjointRASModel, adjointLaminar, dictionary);
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49adjointLaminar::adjointLaminar
50(
51 incompressibleVars& primalVars,
53 objectiveManager& objManager,
54 const word& adjointTurbulenceModelName,
55 const word& modelName
56)
57:
58 adjointRASModel
59 (
60 modelName,
61 primalVars,
62 adjointVars,
63 objManager,
64 adjointTurbulenceModelName
65 )
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
72{
73 const volVectorField& Ua = adjointVars_.Ua();
74 return devReff(Ua);
75}
76
77
79(
80 const volVectorField& U
81) const
82{
84 (
85 "devRhoReff",
95{
96 return
107 (
108 "adjointMeanFlowSource",
110 mesh_,
111 dimensionedVector(dimensionSet(0, 1, -2, 0, 0), Zero)
112 );
113}
114
115
117{
118 // zero contribution
120}
121
126}
127
138 (
139 "adjointEikonalSource" + type(),
150 (
151 "volumeSensTerm" + type(),
153 mesh_,
154 dimensionedTensor(dimensionSet(0, 2, -3, 0, 0), Zero)
155 );
156}
157
158
160(
161 const word& designVarsName
162) const
163{
164 return tmp<scalarField>::New(mesh_.nCells(), Zero);
165}
166
169{
170 // Does nothing. No fields to nullify
171}
172
175{
176 return adjointRASModel::read();
177}
178
179
181{
183}
184
185// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187} // End namespace adjointRASModels
188} // End namespace incompressibleAdjoint
189} // End namespace Foam
190
191// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Manages the adjoint mean flow fields and their mean values.
Abstract base class for incompressible turbulence models.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
virtual bool read()
Read adjointRASProperties dictionary.
Dummy turbulence model for a laminar incompressible flow. Can also be used when the "frozen turbulenc...
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the diffusion term for the momentum equation.
virtual tmp< volTensorField > FISensitivityTerm()
Returns zero field.
virtual void correct()
Correct the primal viscosity field. Redundant?
virtual const boundaryVectorField & adjointMomentumBCSource() const
Returns zero field.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor, i.e. the adjointLaminar stress.
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
virtual const boundaryVectorField & wallShapeSensitivities()
Returns zero field.
virtual tmp< volVectorField > adjointMeanFlowSource()
Source terms to the adjoint momentum equation due to the differentiation of the turbulence model.
virtual const boundaryVectorField & wallFloCoSensitivities()
Returns zero field.
virtual tmp< volScalarField > distanceSensitivities()
Returns zero field.
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const
Returns zero field.
virtual bool read()
Read adjointRASProperties dictionary.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void correct()=0
Solve the adjoint turbulence equations.
Base class for solution control classes.
Class for managing objective functions.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & T
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volVectorField::Boundary boundaryVectorField
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.