objectiveNutSqr.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 
29 #include "objectiveNutSqr.H"
30 #include "createZeroField.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 namespace objectives
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(objectiveNutSqr, 1);
45 (
46  objectiveIncompressible,
47  objectiveNutSqr,
48  dictionary
49 );
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvMesh& mesh,
57  const dictionary& dict,
58  const word& adjointSolverName,
59  const word& primalSolverName
60 )
61 :
62  objectiveIncompressible(mesh, dict, adjointSolverName, primalSolverName),
63  zones_
64  (
65  mesh_.cellZones().indices(this->dict().get<wordRes>("zones"))
66  )
67 {
68  // Allocate source term for the adjoint turbulence model
69  dJdTMvar1Ptr_.reset
70  (
71  createZeroFieldPtr<scalar>
72  (
73  mesh_,
74  "ATMSource",
76  )
77  );
78  // Allocate term to be added to volume-based sensitivity derivatives
79  divDxDbMultPtr_.reset
80  (
81  createZeroFieldPtr<scalar>
82  (
83  mesh_,
84  ("divDxdbMult"+type()) ,
85  //variable dimensions!!
86  //Dummy zeroes. Only the internalField will be used
87  dimless
88  )
89  );
90  // set file pointer
91  //objFunctionFilePtr_ = new OFstream(objFunctionFolder_/type());
92 }
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  J_ = Zero;
100 
102  turbVars = vars_.RASModelVariables();
103  const volScalarField& nut = turbVars->nutRefInst();
104 
105  //scalar zoneVol(Zero);
106  for (const label zI : zones_)
107  {
108  const cellZone& zoneI = mesh_.cellZones()[zI];
109  for (const label cellI : zoneI)
110  {
111  J_ += sqr(nut[cellI])*(mesh_.V()[cellI]);
112  //zoneVol += mesh_.V()[cellI];
113  }
114  }
115  reduce(J_, sumOp<scalar>());
116  //reduce(zoneVol, sumOp<scalar>());
117 
118  return J_;
119 }
120 
121 
123 {
125  turbVars = vars_.RASModelVariables();
126  const singlePhaseTransportModel& lamTransp = vars_.laminarTransport();
127  const volScalarField& nut = turbVars->nutRef();
128 
129  tmp<volScalarField> tnutJacobian = turbVars->nutJacobianVar1(lamTransp);
130  const volScalarField& nutJacobian = tnutJacobian();
131 
133 
134  for (const label zI : zones_)
135  {
136  const cellZone& zoneI = mesh_.cellZones()[zI];
137  for (const label cellI : zoneI)
138  {
139  dJdTMvar1[cellI] = 2.*nut[cellI]*nutJacobian[cellI];
140  }
141  }
142 }
143 
144 
146 {
148  turbVars = vars_.RASModelVariables();
149  const volScalarField& nut = turbVars->nutRef();
150 
151  volScalarField& divDxDbMult = divDxDbMultPtr_();
152 
153  for (const label zI : zones_)
154  {
155  const cellZone& zoneI = mesh_.cellZones()[zI];
156  for (const label cellI : zoneI)
157  {
158  divDxDbMult[cellI] = sqr(nut[cellI]);
159  }
160  }
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace objectives
167 } // End namespace Foam
168 
169 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::objectiveIncompressible::vars_
const incompressibleVars & vars_
Definition: objectiveIncompressible.H:62
Foam::objectives::addToRunTimeSelectionTable
addToRunTimeSelectionTable(objectiveIncompressible, objectiveForce, dictionary)
Foam::objectives::objectiveNutSqr::update_divDxDbMultiplier
void update_divDxDbMultiplier()
Definition: objectiveNutSqr.C:145
Foam::objectives::objectiveNutSqr::J
scalar J()
Return the objective function value.
Definition: objectiveNutSqr.C:97
Foam::objectiveIncompressible::dJdTMvar1
const volScalarField & dJdTMvar1()
Contribution to field adjoint turbulence model variable 1.
Definition: objectiveIncompressible.C:249
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:57
objectiveNutSqr.H
Foam::sumOp
Definition: ops.H:213
Foam::incompressibleVars::RASModelVariables
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Definition: incompressibleVars.C:444
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::objectives::defineTypeNameAndDebug
defineTypeNameAndDebug(objectiveForce, 0)
Foam::objective::J_
scalar J_
Objective function value and weight.
Definition: objective.H:77
createZeroField.H
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::objectives::objectiveNutSqr::update_dJdTMvar1
void update_dJdTMvar1()
Update field to be added to the adjoint turbulence model PDE.
Definition: objectiveNutSqr.C:122
Foam::incompressibleVars::laminarTransport
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
Definition: incompressibleVars.C:418
Foam::objectiveIncompressible
Abstract base class for objective functions in incompressible flows.
Definition: objectiveIncompressible.H:54
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::objectiveIncompressible::dJdTMvar1Ptr_
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
Definition: objectiveIncompressible.H:72
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objectives::objectiveNutSqr::objectiveNutSqr
objectiveNutSqr(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
From components.
Definition: objectiveNutSqr.C:55
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::objective::mesh_
const fvMesh & mesh_
Definition: objective.H:67
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179