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-------------------------------------------------------------------------------
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
32#include "fvPatchFieldMapper.H"
33#include "volFields.H"
34#include "surfaceFields.H"
35#include "emptyFvPatch.H"
36#include "ATCUaGradU.H"
37
38// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39
40Foam::adjointOutletPressureFvPatchScalarField::
41adjointOutletPressureFvPatchScalarField
42(
43 const fvPatch& p,
45)
46:
47 fixedValueFvPatchScalarField(p, iF),
49{}
50
51
52Foam::adjointOutletPressureFvPatchScalarField::
53adjointOutletPressureFvPatchScalarField
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
66Foam::adjointOutletPressureFvPatchScalarField::
67adjointOutletPressureFvPatchScalarField
68(
69 const fvPatch& p,
71 const dictionary& dict
72)
73:
74 fixedValueFvPatchScalarField(p, iF),
75 adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
76{
78 (
79 scalarField("value", dict, p.size())
80 );
81}
82
83
84Foam::adjointOutletPressureFvPatchScalarField::
85adjointOutletPressureFvPatchScalarField
86(
89)
90:
91 fixedValueFvPatchScalarField(tppsf, iF),
93{}
94
95
96// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97
99{
100 if (updated())
101 {
102 return;
103 }
104
105 // Patch normal and surface
106 const scalarField& magSf = patch().magSf();
107 const vectorField nf(patch().nf());
108
109 // Primal flux
110 const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
111 const fvPatchField<vector>& Up = boundaryContrPtr_->Ub();
112
113 // Adjoint velocity
114 const fvPatchField<vector>& Uap = boundaryContrPtr_->Uab();
115 scalarField snGradUan(Uap.snGrad() & nf);
116
117 // Patch normal adjoint velocity
118 scalarField Uap_n(Uap & nf);
119
120 // Patch normal velocity
121 scalarField phiOverSurf(phip/magSf);
122
123 // Momentum diffusion coefficient
124 tmp<scalarField> tmomentumDiffusion =
125 boundaryContrPtr_->momentumDiffusion();
126 const scalarField& momentumDiffusion = tmomentumDiffusion();
127
128 // Part of the diffusive flux related to div(nuEff*dev(grad(Ua).T()))
129 const word& UaName = boundaryContrPtr_->Uab().internalField().name();
130 tmp<tensorField> tgradUab = computePatchGrad<vector>(UaName);
131 const tensorField& gradUab = tgradUab();
132 vectorField explDiffusiveFlux
133 (
134 momentumDiffusion*(gradUab - sphericalTensor::oneThirdI*tr(gradUab))
135 & nf
136 );
137 scalarField normalExplDifFlux(explDiffusiveFlux & nf);
138
139 // Objective function and other explicit contributions
140 tmp<scalarField> tsource = boundaryContrPtr_->pressureSource();
141 scalarField& source = tsource.ref();
142
143 // Contribution from the ATC part (if UaGradU)
144 if (addATCUaGradUTerm())
145 {
146 source += Uap & Up;
147 }
148
149 operator==
150 (
151 (Uap_n*phiOverSurf)
152 + momentumDiffusion*snGradUan
153 + normalExplDifFlux
154 + source
155 );
156
157 fixedValueFvPatchScalarField::updateCoeffs();
158}
159
160
162{
164 writeEntry("value", os);
165 os.writeEntry("solverName", adjointSolverName_);
166}
167
168
169// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170
171namespace Foam
172{
174 (
177 );
178}
179
180// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
static const SphericalTensor oneThirdI
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Base class for solution control classes.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict
Foam::surfaceFields.