ATCUaGradU.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 
30 #include "ATCUaGradU.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(ATCUaGradU, 0);
42 (
43  ATCModel,
44  ATCUaGradU,
45  dictionary
46 );
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 ATCUaGradU::ATCUaGradU
52 (
53  const fvMesh& mesh,
54  const incompressibleVars& primalVars,
55  const incompressibleAdjointVars& adjointVars,
56  const dictionary& dict
57 )
58 :
59  ATCModel(mesh, primalVars, adjointVars, dict)
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 {
67  const volVectorField& U = primalVars_.U();
68  const volVectorField& Ua = adjointVars_.UaInst();
71 
72  // Build Ua to go into the ATC term, based on whether to smooth field or not
73  autoPtr<volVectorField> UaForATC(nullptr);
75  {
76  UaForATC.reset(new volVectorField(fvc::reconstruct(phia)));
77  }
78  else
79  {
80  UaForATC.reset(new volVectorField(Ua));
81  }
82 
83  // Main ATC term
84  ATC_ = -fvc::grad(UaForATC(), "gradUaATC") & U;
85 
86  if (extraConvection_ > 0)
87  {
88  // Implicit part added to increase diagonal dominance
89  // Note: Maybe this needs to be multiplied with the ATClimiter ??
90  UaEqn += extraConvection_*fvm::div(-phi, Ua);
91 
92  // correct rhs due to implicitly augmenting the adjoint convection
93  ATC_ += extraConvection_*(fvc::grad(UaForATC(), "gradUaATC")().T() & U);
94  }
95 
96  // Zero ATC on cells next to given patch types
97  smoothATC();
98 
99  // Actual ATC term
100  UaEqn += fvm::Su(ATC_, Ua);
101 }
102 
103 
105 {
106  tmp<volTensorField> tvolSDTerm
107  (
108  new volTensorField
109  (
110  IOobject
111  (
112  "ATCFISensitivityTerm" + type(),
113  mesh_.time().timeName(),
114  mesh_,
117  ),
118  mesh_,
120  )
121  );
122  volTensorField& volSDTerm = tvolSDTerm.ref();
123 
124  const volVectorField& U = primalVars_.U();
125  const volVectorField& Ua = adjointVars_.Ua();
126 
127  //const volScalarField& mask = getLimiter();
128 
129  volSDTerm -=
130  Ua.component(0) * fvc::grad(U.component(0) * U)
131  + Ua.component(1) * fvc::grad(U.component(1) * U)
132  + Ua.component(2) * fvc::grad(U.component(2) * U);
133 
134  return tvolSDTerm;
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace Foam
141 
142 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:56
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::incompressibleAdjointMeanFlowVars::phiaInst
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
Definition: incompressibleAdjointMeanFlowVars.C:249
Foam::incompressibleVars::phi
const surfaceScalarField & phi() const
Return const reference to volume flux.
Definition: incompressibleVars.C:357
Foam::fvm::Su
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::ATCModel
Base class for selecting the adjoint transpose convection model. Inherits from regIOobject to add loo...
Definition: ATCModel.H:60
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
ATCUaGradU.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::incompressibleAdjointMeanFlowVars::Ua
const volVectorField & Ua() const
Return const reference to velocity.
Definition: incompressibleAdjointMeanFlowVars.C:173
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
Foam::ATCUaGradU::addATC
virtual void addATC(fvVectorMatrix &UaEqn)
Add ATC.
Definition: ATCUaGradU.C:65
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::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::ATCModel::mesh_
const fvMesh & mesh_
Definition: ATCModel.H:79
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::ATCModel::adjointVars_
const incompressibleAdjointVars & adjointVars_
Definition: ATCModel.H:81
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::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::ATCModel::smoothATC
void smoothATC()
Limit ATC field using ATClimiter_.
Definition: ATCModel.C:51
Foam::ATCUaGradU::getFISensitivityTerm
virtual tmp< volTensorField > getFISensitivityTerm() const
Get the FI sensitivity derivatives term coming from the ATC.
Definition: ATCUaGradU.C:104
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::ATCModel::ATC_
volVectorField ATC_
Definition: ATCModel.H:92
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::dimensionedTensor
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedTensor.H:52
Foam::incompressibleAdjointMeanFlowVars::UaInst
const volVectorField & UaInst() const
Return const reference to velocity.
Definition: incompressibleAdjointMeanFlowVars.C:237
Foam::ATCModel::reconstructGradients_
const bool reconstructGradients_
Definition: ATCModel.H:88
Foam::ATCModel::extraConvection_
const scalar extraConvection_
Definition: ATCModel.H:85
Foam::ATCModel::primalVars_
const incompressibleVars & primalVars_
Definition: ATCModel.H:80
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54