DAObjFuncMass.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAObjFuncMass.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAObjFuncMass, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncMass, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const fvMesh& mesh,
21  const DAOption& daOption,
22  const DAModel& daModel,
23  const DAIndex& daIndex,
24  const DAResidual& daResidual,
25  const word objFuncName,
26  const word objFuncPart,
27  const dictionary& objFuncDict)
28  : DAObjFunc(
29  mesh,
30  daOption,
31  daModel,
32  daIndex,
33  daResidual,
34  objFuncName,
35  objFuncPart,
36  objFuncDict)
37 {
38  // Assign type, this is common for all objectives
39  objFuncDict_.readEntry<word>("type", objFuncType_);
40 
41  // setup the connectivity for mass. It does actually not depends on D
42  // here we set zero level of D as a place holder.
43  objFuncConInfo_ = {{"D"}};
44 
45  objFuncDict_.readEntry<scalar>("scale", scale_);
46 }
47 
50  const labelList& objFuncFaceSources,
51  const labelList& objFuncCellSources,
52  scalarList& objFuncFaceValues,
53  scalarList& objFuncCellValues,
54  scalar& objFuncValue)
55 {
56  /*
57  Description:
58  Calculate the mass = density * volume of the mesh
59 
60  Input:
61  objFuncFaceSources: List of face source (index) for this objective
62 
63  objFuncCellSources: List of cell source (index) for this objective
64 
65  Output:
66  objFuncFaceValues: the discrete value of objective for each face source (index).
67  This will be used for computing df/dw in the adjoint.
68 
69  objFuncCellValues: the discrete value of objective on each cell source (index).
70  This will be used for computing df/dw in the adjoint.
71 
72  objFuncValue: the sum of objective, reduced across all processors and scaled by "scale"
73  */
74 
75  // initialize faceValues to zero
76  forAll(objFuncCellValues, idxI)
77  {
78  objFuncCellValues[idxI] = 0.0;
79  }
80  // initialize objFunValue
81  objFuncValue = 0.0;
82 
83  const objectRegistry& db = mesh_.thisDb();
84  const volScalarField& rho = db.lookupObject<volScalarField>("solid:rho");
85 
86  // calculate mass
87  forAll(objFuncCellSources, idxI)
88  {
89  const label& cellI = objFuncCellSources[idxI];
90  scalar volume = mesh_.V()[cellI];
91  objFuncCellValues[idxI] = scale_ * volume * rho[cellI];
92  objFuncValue += objFuncCellValues[idxI];
93  }
94 
95  // need to reduce the sum of force across all processors
96  reduce(objFuncValue, sumOp<scalar>());
97 
98  return;
99 }
100 
101 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
102 
103 } // End namespace Foam
104 
105 // ************************************************************************* //
Foam::DAObjFunc::objFuncType_
word objFuncType_
the type of the objective function
Definition: DAObjFunc.H:66
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DAOption
Definition: DAOption.H:29
daOption
DAOption daOption(mesh, pyOptions_)
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAObjFunc::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAObjFunc.H:45
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAObjFuncMass::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncMass.C:49
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAObjFuncMass::DAObjFuncMass
DAObjFuncMass(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAResidual &daResidual, const word objFuncName, const word objFuncPart, const dictionary &objFuncDict)
Definition: DAObjFuncMass.C:19
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DAObjFunc::scale_
scalar scale_
scale of the objective function
Definition: DAObjFunc.H:87
daModel
DAModel daModel(mesh, daOption)
Foam::DAObjFunc::objFuncConInfo_
List< List< word > > objFuncConInfo_
the connectivity information for the objective function
Definition: DAObjFunc.H:93
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAObjFunc::objFuncDict_
const dictionary & objFuncDict_
dictionary containing the information for the objective function
Definition: DAObjFunc.H:69
Foam::DAObjFunc
Definition: DAObjFunc.H:33
DAObjFuncMass.H
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8