DAIndex.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  Handel state indexing in serial and in parallel
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DAIndex_H
12 #define DAIndex_H
13 
14 #include "fvOptions.H"
15 #include "surfaceFields.H"
16 #include "DAOption.H"
17 #include "DAUtility.H"
18 #include "DAStateInfo.H"
19 #include "DAModel.H"
20 #include "globalIndex.H"
21 #include "DAMacroFunctions.H"
22 
23 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
24 
25 namespace Foam
26 {
27 
28 /*---------------------------------------------------------------------------*\
29  Class DAIndex Declaration
30 \*---------------------------------------------------------------------------*/
31 
32 class DAIndex
33 {
34 
35 private:
37  DAIndex(const DAIndex&);
38 
40  void operator=(const DAIndex&);
41 
42 protected:
44  const fvMesh& mesh_;
45 
48 
50  const DAModel& daModel_;
51 
53  HashTable<wordList> stateInfo_;
54 
56  void writeAdjointIndexing();
57 
58 public:
60  DAIndex(
61  const fvMesh& mesh,
62  const DAOption& daOption,
63  const DAModel& daModel);
64 
66  virtual ~DAIndex()
67  {
68  }
69 
70  // Members
71 
72  // adjoint indexing info
74  List<word> adjStateNames;
75 
77  HashTable<word> adjStateType;
78 
79  // mesh sizes
80  // local sizes
81 
83  label nLocalCells;
84 
86  label nLocalFaces;
87 
89  label nLocalPoints;
90 
92  label nLocalXv;
93 
96 
99 
102 
105 
110  labelList bFacePatchI;
111 
116  labelList bFaceFaceI;
117 
122  HashTable<label> stateLocalIndexOffset;
123 
125  HashTable<label> adjStateID;
126 
127  // glocal sizes
130 
133 
135  label nGlobalXv;
136 
139 
142 
143  // adjoint sizes
146 
149 
152 
155 
158 
161 
164 
166 
168  globalIndex globalCellNumbering;
169  globalIndex globalCellVectorNumbering; // similar to globalCellNumbering but has 3 components per cell
170  globalIndex globalFaceNumbering;
173  globalIndex globalXvNumbering;
178 
180  labelIOList pointProcAddressing;
181 
184 
186  labelList isCoupledFace;
187 
190 
193 
195  labelList faceOwner;
196 
199 
201  labelList phiLocalOffset;
202 
203  // Member functions
204 
206  void calcStateLocalIndexOffset(HashTable<label>& offset);
207 
209  void calcAdjStateID(HashTable<label>& adjStateID);
210 
212  void calcLocalIdxLists(
213  wordList& adjStateName4LocalAdjIdx,
214  scalarList& cellIFaceI4LocalAdjIdx);
215 
218  const word stateName,
219  const label idxI,
220  const label comp = -1) const;
221 
224  const word stateName,
225  const label idxI,
226  const label comp = -1) const;
227 
229  label getGlobalXvIndex(
230  const label idxPoint,
231  const label idxCoord) const;
232 
234  label getLocalXvIndex(
235  const label idxPoint,
236  const label idxCoord) const;
237 
239  label getLocalCellIndex(const label cellI) const;
240 
242  label getGlobalCellIndex(const label cellI) const;
243 
246  const label cellI,
247  const label comp) const;
248 
251  const label cellI,
252  const label comp) const;
253 
255  label getLocalFaceIndex(const label faceI) const;
256 
258  label getGlobalFaceIndex(const label faceI) const;
259 
261  void calcAdjStateID4GlobalAdjIdx(labelList& adjStateID4GlobalAdjIdx) const;
262 
264  void printIndices() const
265  {
266  //Info<<"adjStateNames"<<adjStateNames<<endl;
267  Info << "Adjoint States: " << adjStateType << endl;
268 
269  // Print relevant sizes to screen
270  Info << "Global Cells: " << nGlobalCells << endl;
271  Info << "Global Faces: " << nGlobalFaces << endl;
272  Info << "Global Xv: " << nGlobalXv << endl;
273  Info << "Undecomposed points: " << nUndecomposedPoints << endl;
274  Info << "Global Adjoint States: " << nGlobalAdjointStates << endl;
275  }
276 
278  void printMatChars(const Mat matIn) const;
279 
281  void getMatNonZeros(
282  const Mat matIn,
283  label& maxCols,
284  scalar& allNonZeros) const;
285 };
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #endif
294 
295 // ************************************************************************* //
Foam::DAIndex::getGlobalFaceIndex
label getGlobalFaceIndex(const label faceI) const
get global face index for a local face index
Definition: DAIndex.C:886
Foam::DAIndex::faceOwner
labelList faceOwner
owner cell of a given face
Definition: DAIndex.H:195
DAOption.H
Foam::DAIndex::globalFaceNumbering
globalIndex globalFaceNumbering
Definition: DAIndex.H:170
Foam::DAIndex::~DAIndex
virtual ~DAIndex()
Destructor.
Definition: DAIndex.H:66
Foam::DAIndex::globalVolScalarFieldNumbering
globalIndex globalVolScalarFieldNumbering
Definition: DAIndex.H:174
Foam::DAIndex::nLocalFaces
label nLocalFaces
local face size
Definition: DAIndex.H:86
Foam::DAIndex::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAIndex.H:50
Foam::DAIndex::getMatNonZeros
void getMatNonZeros(const Mat matIn, label &maxCols, scalar &allNonZeros) const
get number of non zeros for a matrix
Definition: DAIndex.C:1236
Foam::DAIndex::globalAdjointStateNumbering
globalIndex globalAdjointStateNumbering
Definition: DAIndex.H:167
Foam::DAIndex::getLocalCellIndex
label getLocalCellIndex(const label cellI) const
get local cell index for a local cell index
Definition: DAIndex.C:760
Foam::DAIndex::nLocalXv
label nLocalXv
local Xv size (point size*3)
Definition: DAIndex.H:92
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAIndex::nUndecomposedPoints
label nUndecomposedPoints
number of points for the un-decomposed domain
Definition: DAIndex.H:183
Foam::DAIndex::adjStateName4LocalAdjIdx
wordList adjStateName4LocalAdjIdx
given a local adjoint state index, return its state name
Definition: DAIndex.H:189
Foam::DAOption
Definition: DAOption.H:29
Foam::DAIndex::getGlobalAdjointStateIndex
label getGlobalAdjointStateIndex(const word stateName, const label idxI, const label comp=-1) const
get global adjoint index for a given state name, cell/face indxI and its component (optional,...
Definition: DAIndex.C:661
DAUtility.H
Foam::DAIndex::nGlobalXv
label nGlobalXv
global Xv size (global face size*3)
Definition: DAIndex.H:135
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAIndex::nLocalAdjointBoundaryStates
label nLocalAdjointBoundaryStates
local adjoint boundary states that do not contain phis
Definition: DAIndex.H:151
Foam::DAIndex::cellIFaceI4LocalAdjIdx
scalarList cellIFaceI4LocalAdjIdx
given a local adjoint state index, return its cell/face index
Definition: DAIndex.H:192
Foam::DAIndex::nGlobalObjFuncFaces
label nGlobalObjFuncFaces
global objective function face size
Definition: DAIndex.H:138
Foam::DAIndex::globalCellVectorNumbering
globalIndex globalCellVectorNumbering
Definition: DAIndex.H:169
Foam::DAIndex::adjStateType
HashTable< word > adjStateType
hash table of adjoint state types, e.g., volVectorState for a given state name
Definition: DAIndex.H:77
DAMacroFunctions.H
Foam::DAIndex::nLocalCells
label nLocalCells
local cell size
Definition: DAIndex.H:83
Foam::DAIndex::isCoupledFace
labelList isCoupledFace
for a given face index, return whether this face is a coupled boundary face
Definition: DAIndex.H:186
Foam::DAIndex::writeAdjointIndexing
void writeAdjointIndexing()
write the adjoint indexing for debugging
Definition: DAIndex.C:1283
Foam::DAIndex::phiLocalOffset
labelList phiLocalOffset
phi local indexing offset for cell-by-cell indexing
Definition: DAIndex.H:201
Foam::DAIndex::stateInfo_
HashTable< wordList > stateInfo_
the StateInfo_ list from DAStateInfo object
Definition: DAIndex.H:53
Foam::DAIndex::calcLocalIdxLists
void calcLocalIdxLists(wordList &adjStateName4LocalAdjIdx, scalarList &cellIFaceI4LocalAdjIdx)
compute local lists such as adjStateName4LocalAdjIdx and cellIFaceI4LocalAdjIdx;
Definition: DAIndex.C:441
Foam::DAIndex::calcStateLocalIndexOffset
void calcStateLocalIndexOffset(HashTable< label > &offset)
calculate stateLocalIndexOffset
Definition: DAIndex.C:186
Foam::DAIndex::nLocalBoundaryFaces
label nLocalBoundaryFaces
local boundary face size
Definition: DAIndex.H:95
Foam::DAIndex::pointProcAddressing
labelIOList pointProcAddressing
a list to map the point index from decomposed domains to the original un-decomposed domain
Definition: DAIndex.H:180
DAModel.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex::nVolVectorStates
label nVolVectorStates
number of state variables for volVectorField
Definition: DAIndex.H:157
Foam::DAIndex::globalXvNumbering
globalIndex globalXvNumbering
Definition: DAIndex.H:173
Foam::DAIndex::globalCellNumbering
globalIndex globalCellNumbering
Definition: DAIndex.H:168
Foam::DAIndex::globalCoupledBFaceNumbering
globalIndex globalCoupledBFaceNumbering
Definition: DAIndex.H:171
Foam::DAIndex::nVolScalarStates
label nVolScalarStates
number of state variables for volScalarField
Definition: DAIndex.H:154
Foam::DAIndex::getGlobalCellVectorIndex
label getGlobalCellVectorIndex(const label cellI, const label comp) const
get global cell index for a local cell vector index
Definition: DAIndex.C:834
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAIndex::adjStateNames
List< word > adjStateNames
list of adjoint state names for a specific solver
Definition: DAIndex.H:74
Foam::DAIndex::nLocalBoundaryPatches
label nLocalBoundaryPatches
local boundary patch size
Definition: DAIndex.H:101
Foam::DAIndex::getGlobalXvIndex
label getGlobalXvIndex(const label idxPoint, const label idxCoord) const
get global Xv index for a given point index and coordinate component (x, y, or z)
Definition: DAIndex.C:702
Foam::DAIndex::phiAccumulatdOffset
labelList phiAccumulatdOffset
the accumulated phi indexing offset for cell-by-cell indexing
Definition: DAIndex.H:198
Foam::DAIndex::nGlobalCells
label nGlobalCells
global cell size
Definition: DAIndex.H:129
Foam::DAModel
Definition: DAModel.H:59
Foam::DAIndex::globalVolVectorFieldNumbering
globalIndex globalVolVectorFieldNumbering
Definition: DAIndex.H:175
Foam::DAIndex::calcAdjStateID
void calcAdjStateID(HashTable< label > &adjStateID)
set adjoint state unique ID: adjStateID
Definition: DAIndex.C:397
Foam::DAIndex::stateLocalIndexOffset
HashTable< label > stateLocalIndexOffset
Definition: DAIndex.H:122
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAIndex::globalSurfaceScalarFieldNumbering
globalIndex globalSurfaceScalarFieldNumbering
Definition: DAIndex.H:176
Foam::DAIndex::getLocalFaceIndex
label getLocalFaceIndex(const label faceI) const
get local face index for a local face index
Definition: DAIndex.C:864
Foam::DAIndex::nGlobalCoupledBFaces
label nGlobalCoupledBFaces
global coupled boundary face size
Definition: DAIndex.H:141
Foam::DAIndex::nLocalPoints
label nLocalPoints
local point size
Definition: DAIndex.H:89
Foam::DAIndex::adjStateID
HashTable< label > adjStateID
a unique number ID for adjoint states, it depends on the sequence of adjStateNames
Definition: DAIndex.H:125
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAIndex::getLocalXvIndex
label getLocalXvIndex(const label idxPoint, const label idxCoord) const
get local Xv index for a given point index and coordinate component (x, y, or z)
Definition: DAIndex.C:733
Foam::DAIndex::calcAdjStateID4GlobalAdjIdx
void calcAdjStateID4GlobalAdjIdx(labelList &adjStateID4GlobalAdjIdx) const
compute global list adjStateID4GlobalAdjIdx that stores the stateID for a given globalAdjIndx
Definition: DAIndex.C:912
Foam::DAIndex::getLocalCellVectorIndex
label getLocalCellVectorIndex(const label cellI, const label comp) const
get local cell index for a local cell vector index
Definition: DAIndex.C:808
Foam::DAIndex::printIndices
void printIndices() const
print all the sizes and indices
Definition: DAIndex.H:264
Foam::DAIndex::getGlobalCellIndex
label getGlobalCellIndex(const label cellI) const
get global cell index for a local cell index
Definition: DAIndex.C:782
DAStateInfo.H
Foam::DAIndex::mesh_
const fvMesh & mesh_
Foam::fvMesh object.
Definition: DAIndex.H:44
daModel
DAModel daModel(mesh, daOption)
Foam::DAIndex::nGlobalFaces
label nGlobalFaces
global face size
Definition: DAIndex.H:132
Foam::DAIndex::printMatChars
void printMatChars(const Mat matIn) const
print petsc matrix statistics such as the max on/off diagonral ratio
Definition: DAIndex.C:1010
Foam::DAIndex::getLocalAdjointStateIndex
label getLocalAdjointStateIndex(const word stateName, const label idxI, const label comp=-1) const
get local adjoint index for a given state name, cell/face indxI and its component (optional,...
Definition: DAIndex.C:516
Foam::DAIndex::nModelStates
label nModelStates
number of model states, NOTE: they are counted separately
Definition: DAIndex.H:163
Foam::DAIndex::nGlobalAdjointStates
label nGlobalAdjointStates
number of global adjoint states (including all cells and faces)
Definition: DAIndex.H:148
Foam::DAIndex::nLocalCoupledBFaces
label nLocalCoupledBFaces
local coupled boundary patch size
Definition: DAIndex.H:104
Foam::DAIndex::globalObjFuncGeoNumbering
globalIndex globalObjFuncGeoNumbering
Definition: DAIndex.H:172
Foam::DAIndex::daOption_
const DAOption & daOption_
Foam::DAOption object.
Definition: DAIndex.H:47
Foam::DAIndex::nSurfaceScalarStates
label nSurfaceScalarStates
number of state variables for surfaceScalarField
Definition: DAIndex.H:160
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116
Foam::DAIndex::nLocalAdjointStates
label nLocalAdjointStates
number of local adjoint states (including all cells and faces)
Definition: DAIndex.H:145