adjointOutletKaFvPatchScalarField.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-2020 PCOpt/NTUA
9  Copyright (C) 2013-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 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedValueFvPatchScalarField(p, iF),
49 {}
50 
51 
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61  adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
62 {}
63 
64 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  fixedValueFvPatchScalarField(p, iF),
74 {
75  this->readValueEntry(dict, IOobjectOption::MUST_READ);
76 }
77 
78 
80 (
83 )
84 :
85  fixedValueFvPatchScalarField(tppsf, iF),
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  if (updated())
95  {
96  return;
97  }
98  vectorField nf(patch().nf());
99 
100  const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
101  tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
102  const scalarField& nuEff = tnuEff();
103 
104  // Patch-adjacent ka
105  tmp<scalarField> tkaNei(patchInternalField());
106  const scalarField& kaNei = tkaNei();
107 
108  const scalarField& delta = patch().deltaCoeffs();
109 
110  // Source from the objective
111  tmp<scalarField> source = boundaryContrPtr_->adjointTMVariable1Source();
112 
113  operator==
114  (
115  (nuEff*delta*kaNei - source)
116  /((Ub & nf) + nuEff*delta)
117  );
118 
119  fixedValueFvPatchScalarField::updateCoeffs();
120 }
121 
122 
124 {
126  os.writeEntry("solverName", adjointSolverName_);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 makePatchTypeField(fvPatchScalarField, adjointOutletKaFvPatchScalarField);
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace Foam
138 
139 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
dictionary dict
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
autoPtr< boundaryAdjointContribution > boundaryContrPtr_
Engine to manage contributions of the objective functions to the adjoint boundary conditions...
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
adjointOutletKaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Base class for solution control classes.
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.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
word adjointSolverName_
adjointSolver name corresponding to field
Namespace for OpenFOAM.