adjointFarFieldNuaTildaFvPatchScalarField.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 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
48 )
49 :
50  fixedValueFvPatchScalarField(p, iF),
52 {}
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
65  adjointScalarBoundaryCondition(p, iF, ptf.adjointSolverName_)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  fixedValueFvPatchScalarField(p, iF),
78  adjointScalarBoundaryCondition(p, iF, dict.get<word>("solverName"))
79 {
80  this->readValueEntry(dict, IOobjectOption::MUST_READ);
81 }
82 
83 
86 (
89 )
90 :
91  fixedValueFvPatchScalarField(tppsf, iF),
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (updated())
101  {
102  return;
103  }
104  vectorField nf(patch().nf());
105 
106  const fvPatchField<vector>& Ub = boundaryContrPtr_->Ub();
107  tmp<scalarField> tnuEff(boundaryContrPtr_->TMVariable1Diffusion());
108  const scalarField& nuEff = tnuEff();
109  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
110 
111  // Patch-adjacent nuaTilda nuaTildaNei
112  tmp<scalarField> tnuaTildaNei(patchInternalField());
113  const scalarField& nuaTildaNei = tnuaTildaNei();
114 
115  // Patch deltas
116  const scalarField& delta = patch().deltaCoeffs();
117 
118  operator==
119  (
120  pos(phip)
121  *(
122  (nuEff*delta*nuaTildaNei)
123  /((Ub & nf) + nuEff*delta)
124  )
125  );
127  fixedValueFvPatchScalarField::updateCoeffs();
128 }
129 
130 
133 (
134  const tmp<scalarField>&
135 ) const
136 {
137  const fvsPatchField<scalar>& phip = boundaryContrPtr_->phib();
138 
139  // For fixedValue nuTilda patches
141 }
142 
143 
146 (
147  const tmp<scalarField>&
148 ) const
149 {
151 
152  // For zeroGradient nuTilda patches
153  return tmp<Field<scalar>>::New(pos(phip)*(*this));
154 }
155 
156 
158 {
160  os.writeEntry("solverName", adjointSolverName_);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
168 (
170  adjointFarFieldNuaTildaFvPatchScalarField
171 );
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 } // End namespace Foam
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.
scalar delta
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< Field< scalar > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual tmp< Field< scalar > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
dimensionedScalar pos(const dimensionedScalar &ds)
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.
adjointFarFieldNuaTildaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const word null
An empty word.
Definition: word.H:84
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
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.
adjointBoundaryCondition< scalar > adjointScalarBoundaryCondition
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
word adjointSolverName_
adjointSolver name corresponding to field
Namespace for OpenFOAM.