boundaryAdjointContribution.H
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 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 Class
30  Foam::boundaryAdjointContribution
31 
32 Description
33  Abstract base class for computing contributions of the objective functions
34  to the adjoint boundary conditions
35 
36 SourceFiles
37  boundaryAdjointContribution.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef boundaryAdjointContribution_H
42 #define boundaryAdjointContribution_H
43 
44 #include "IOdictionary.H"
45 #include "autoPtr.H"
46 #include "runTimeSelectionTables.H"
47 #include "fvPatchFields.H"
48 #include "fvsPatchFields.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class boundaryAdjointContribution Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 {
61 private:
62 
63  // Private Member Functions
64 
65  //- No copy construct
67  (
69  ) = delete;
70 
71  //- No copy assignment
72  void operator=(const boundaryAdjointContribution&) = delete;
73 
74 
75 protected:
76 
77  // Protected data
78 
79  const fvPatch& patch_;
80 
81 
82  // Protected Member Functions
83 
84  template<class returnType, class sourceType, class castType>
86  (
87  PtrList<sourceType>& sourceList,
88  const fvPatchField<returnType>&(castType::*boundaryFunction)
89  (const label),
90  bool (castType::*hasFunction)() const
91  );
92 
93 
94 public:
95 
96  //- Runtime type information
97  TypeName("boundaryAdjointContribution");
98 
99 
100  // Declare run-time constructor selection table
101 
103  (
104  autoPtr,
106  dictionary,
107  (
108  const word& managerName,
109  const word& adjointSolverName,
110  const word& simulationType,
111  const fvPatch& patch
112  ),
113  (managerName, adjointSolverName, simulationType, patch)
114  );
115 
116  // Constructors
117 
118  //- Construct from components
120  (
121  const word& managerName,
122  const word& adjointSolverName,
123  const word& simulationType,
124  const fvPatch& patch
125  );
126 
127 
128  // Selectors
129 
130  //- Return a reference to the selected turbulence model
132  (
133  const word& managerName,
134  const word& adjointSolverName,
135  const word& simulationType,
136  const fvPatch& patch
137  );
138 
139 
140  //- Destructor
141  virtual ~boundaryAdjointContribution() = default;
142 
143 
144  // Member Functions
145 
146  // Contribution to surface sensitivities for a specific patch
147  virtual tmp<scalarField> pressureSource() = 0;
148  virtual tmp<vectorField> velocitySource() = 0;
153  virtual tmp<scalarField> dJdnut();
154  virtual tmp<tensorField> dJdGradU();
155  virtual tmp<scalarField> energySource() = 0;
156 
157  virtual tmp<scalarField> momentumDiffusion() = 0;
158  virtual tmp<scalarField> laminarDiffusivity() = 0;
159  virtual tmp<scalarField> thermalDiffusion() = 0;
160  virtual tmp<scalarField> wallDistance() = 0;
161 
164  virtual tmp<scalarField> TMVariable1();
165  virtual tmp<scalarField> TMVariable2();
166 
167  // References to primal and adjoint fields for the specific patch
168  virtual const fvPatchVectorField& Ub() const = 0;
169  virtual const fvPatchScalarField& pb() const = 0;
170  virtual const fvsPatchScalarField& phib() const = 0;
172  virtual const fvPatchVectorField& Uab() const = 0;
173  virtual const fvPatchScalarField& pab() const = 0;
174  virtual const fvsPatchScalarField& phiab() const = 0;
175 
176  // Field suffixes for primal and adjoint fields
177  virtual const word primalSolverName() const = 0;
178  virtual const word adjointSolverName() const = 0;
179 };
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace Foam
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #ifdef NoRepository
190 #endif
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 #endif
195 
196 // ************************************************************************* //
virtual tmp< scalarField > momentumDiffusion()=0
virtual const word adjointSolverName() const =0
virtual const fvPatchVectorField & Uab() const =0
virtual tmp< fvPatchScalarField > turbulentDiffusivity() const
virtual tmp< scalarField > pressureSource()=0
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static autoPtr< boundaryAdjointContribution > New(const word &managerName, const word &adjointSolverName, const word &simulationType, const fvPatch &patch)
Return a reference to the selected turbulence model.
virtual tmp< scalarField > energySource()=0
TypeName("boundaryAdjointContribution")
Runtime type information.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual const fvPatchScalarField & pab() const =0
virtual const fvPatchScalarField & pb() const =0
virtual tmp< scalarField > laminarDiffusivity()=0
virtual const fvsPatchScalarField & phiab() const =0
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< scalarField > adjointTMVariable2Source()
declareRunTimeSelectionTable(autoPtr, boundaryAdjointContribution, dictionary,(const word &managerName, const word &adjointSolverName, const word &simulationType, const fvPatch &patch),(managerName, adjointSolverName, simulationType, patch))
virtual tmp< scalarField > adjointTMVariable1Source()
Abstract base class for computing contributions of the objective functions to the adjoint boundary co...
virtual ~boundaryAdjointContribution()=default
Destructor.
virtual const fvPatchVectorField & Ub() const =0
virtual const word primalSolverName() const =0
virtual tmp< scalarField > thermalDiffusion()=0
virtual tmp< scalarField > TMVariable2Diffusion()
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual const fvsPatchScalarField & phib() const =0
const std::string patch
OpenFOAM patch number as a std::string.
virtual tmp< scalarField > wallDistance()=0
virtual tmp< vectorField > velocitySource()=0
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
tmp< Field< returnType > > sumContributions(PtrList< sourceType > &sourceList, const fvPatchField< returnType > &(castType::*boundaryFunction)(const label), bool(castType::*hasFunction)() const)
virtual tmp< vectorField > normalVelocitySource()=0
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual tmp< scalarField > TMVariable1Diffusion()
Namespace for OpenFOAM.
virtual tmp< vectorField > tangentVelocitySource()=0