DAFunctionPatchMean.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAFunctionPatchMean.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAFunctionPatchMean, 0);
16 addToRunTimeSelectionTable(DAFunction, DAFunctionPatchMean, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const fvMesh& mesh,
21  const DAOption& daOption,
22  const DAModel& daModel,
23  const DAIndex& daIndex,
24  const word functionName)
25  : DAFunction(
26  mesh,
27  daOption,
28  daModel,
29  daIndex,
30  functionName)
31 {
32  functionDict_.readEntry<word>("varName", varName_);
33 
34  functionDict_.readEntry<word>("varType", varType_);
35 
36  functionDict_.readEntry<label>("index", index_);
37 }
38 
41 {
42  /*
43  Description:
44  Calculate the patch mean
45  */
46 
47  // calculate the area of all the patches. We need to recompute because the surface area
48  // may change during the optimization
49  areaSum_ = 0.0;
50  forAll(faceSources_, idxI)
51  {
52  const label& functionFaceI = faceSources_[idxI];
53  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
54  const label patchI = daIndex_.bFacePatchI[bFaceI];
55  const label faceI = daIndex_.bFaceFaceI[bFaceI];
56  areaSum_ += mesh_.magSf().boundaryField()[patchI][faceI];
57  }
58  reduce(areaSum_, sumOp<scalar>());
59 
60  // initialize objFunValue
61  scalar functionValue = 0.0;
62 
63  const objectRegistry& db = mesh_.thisDb();
64 
65  if (varType_ == "scalar")
66  {
67  const volScalarField& var = db.lookupObject<volScalarField>(varName_);
68 
69  // calculate area weighted heat flux
70  forAll(faceSources_, idxI)
71  {
72  const label& functionFaceI = faceSources_[idxI];
73  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
74  const label patchI = daIndex_.bFacePatchI[bFaceI];
75  const label faceI = daIndex_.bFaceFaceI[bFaceI];
76  scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
77  functionValue += scale_ * area * var.boundaryField()[patchI][faceI] / areaSum_;
78  }
79  }
80  else if (varType_ == "vector")
81  {
82  const volVectorField& var = db.lookupObject<volVectorField>(varName_);
83 
84  // calculate area weighted heat flux
85  forAll(faceSources_, idxI)
86  {
87  const label& functionFaceI = faceSources_[idxI];
88  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
89  const label patchI = daIndex_.bFacePatchI[bFaceI];
90  const label faceI = daIndex_.bFaceFaceI[bFaceI];
91  scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
92  functionValue += scale_ * area * var.boundaryField()[patchI][faceI][index_] / areaSum_;
93  }
94  }
95  else
96  {
97  FatalErrorIn("DAFunctionPatchMean::calcFunction")
98  << "varType not valid. Options are scalar or vector"
99  << abort(FatalError);
100  }
101 
102  // need to reduce the sum of force across all processors
103  reduce(functionValue, sumOp<scalar>());
104 
105  // check if we need to calculate refDiff.
106  this->calcRefVar(functionValue);
107 
108  return functionValue;
109 }
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 } // End namespace Foam
114 
115 // ************************************************************************* //
Foam::DAFunction::faceSources_
labelList faceSources_
a sorted list of all face sources for the objective function
Definition: DAFunction.H:67
Foam::DAFunctionPatchMean::varType_
word varType_
type of the variable either vector or scalar
Definition: DAFunctionPatchMean.H:38
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAFunction::calcRefVar
void calcRefVar(scalar &functionValue)
calculate (var-ref)^2
Definition: DAFunction.C:216
Foam::DAFunctionPatchMean::varName_
word varName_
name of the variable
Definition: DAFunctionPatchMean.H:35
Foam::DAFunction::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAFunction.H:52
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
DAFunctionPatchMean.H
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:57
Foam::DAFunction
Definition: DAFunction.H:31
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam::DAFunctionPatchMean::calcFunction
virtual scalar calcFunction()
calculate the value of objective function
Definition: DAFunctionPatchMean.C:40
Foam
Definition: checkGeometry.C:32
Foam::DAFunctionPatchMean::DAFunctionPatchMean
DAFunctionPatchMean(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const word functionName)
Definition: DAFunctionPatchMean.C:19
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAFunctionPatchMean::index_
label index_
if vector which element/index?
Definition: DAFunctionPatchMean.H:41
Foam::DAFunction::functionDict_
dictionary functionDict_
dictionary containing the information for the objective function
Definition: DAFunction.H:64
Foam::DAFunction::scale_
scalar scale_
scale of the objective function
Definition: DAFunction.H:73
Foam::DAFunctionPatchMean::areaSum_
scalar areaSum_
the area of all total pressure patches
Definition: DAFunctionPatchMean.H:32
Foam::DAFunction::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAFunction.H:43
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116