adjointOutletPressureFvPatchScalarField.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  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 
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "emptyFvPatch.H"
36 #include "ATCUaGradU.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedValueFvPatchScalarField(p, iF),
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
62  adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
63 {}
64 
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  fixedValueFvPatchScalarField(p, iF),
75  adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
76 {
77  this->readValueEntry(dict, IOobjectOption::MUST_READ);
78 }
79 
80 
83 (
86 )
87 :
88  fixedValueFvPatchScalarField(tppsf, iF),
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (updated())
98  {
99  return;
100  }
101 
102  // Patch normal and surface
103  const scalarField& magSf = patch().magSf();
104  const vectorField nf(patch().nf());
105 
106  // Primal flux
107  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
108  const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
109 
110  // Adjoint velocity
111  const fvPatchField<vector>& Uap = boundaryContrPtr_->Uab();
112  scalarField snGradUan(Uap.snGrad() & nf);
113 
114  // Patch normal adjoint velocity
115  scalarField Uap_n(Uap & nf);
116 
117  // Patch normal velocity
118  scalarField phiOverSurf(phip/magSf);
119 
120  // Momentum diffusion coefficient
121  tmp<scalarField> tmomentumDiffusion =
122  boundaryContrPtr_->momentumDiffusion();
123  const scalarField& momentumDiffusion = tmomentumDiffusion();
124 
125  // Part of the diffusive flux related to div(nuEff*dev(grad(Ua).T()))
126  const word& UaName = boundaryContrPtr_->Uab().internalField().name();
127  tmp<tensorField> tgradUab = computePatchGrad<vector>(UaName);
128  const tensorField& gradUab = tgradUab();
129  vectorField explDiffusiveFlux
130  (
131  momentumDiffusion*(gradUab - sphericalTensor::oneThirdI*tr(gradUab))
132  & nf
133  );
134  scalarField normalExplDifFlux(explDiffusiveFlux & nf);
135 
136  // Objective function and other explicit contributions
137  tmp<scalarField> tsource = boundaryContrPtr_->pressureSource();
138  scalarField& source = tsource.ref();
139 
140  // Contribution from the ATC part (if UaGradU)
141  if (addATCUaGradUTerm())
142  {
143  source += Uap & Up;
144  }
145 
146  operator==
147  (
148  (Uap_n*phiOverSurf)
149  + momentumDiffusion*snGradUan
150  + normalExplDifFlux
151  + source
152  );
153 
154  fixedValueFvPatchScalarField::updateCoeffs();
155 }
156 
157 
159 {
161  os.writeEntry("solverName", adjointSolverName_);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 namespace Foam
169 {
171  (
173  adjointOutletPressureFvPatchScalarField
174  );
175 }
176 
177 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
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
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
adjointOutletPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal 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
static const SphericalTensor oneThirdI
bool addATCUaGradUTerm()
Whether to add the extra term from the UaGradU formulation.
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
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)
Namespace for OpenFOAM.