DAUtility.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  DAUtility contains basic functions such as matrix file IO, and
8  it is independent of fvMesh. All the functions in DAUtility are static,
9  so they can be directly called, e.g., DAUtility::boundVar(...)
10 
11 \*---------------------------------------------------------------------------*/
12 
13 #ifndef DAUtility_H
14 #define DAUtility_H
15 
16 #include <petscksp.h>
17 #include "Python.h"
18 #include "fvOptions.H"
19 #include "globalIndex.H"
20 #include "IOMRFZoneListDF.H"
21 
22 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
23 
24 namespace Foam
25 {
26 
27 typedef void (*pyComputeInterface)(const double*, int, double*, int, void*);
28 typedef void (*pyJacVecProdInterface)(const double*, double*, int, const double*, const double*, int, void*);
29 typedef void (*pySetCharInterface)(const char*, void*);
30 
31 /*---------------------------------------------------------------------------*\
32  Class DAUtility Declaration
33 \*---------------------------------------------------------------------------*/
34 class DAUtility
35 {
36 
37 private:
39  DAUtility(const DAUtility&);
40 
42  void operator=(const DAUtility&);
43 
44 public:
46  DAUtility();
47 
49  virtual ~DAUtility();
50 
52  static void pyDict2OFDict(
53  PyObject* pyDict,
54  dictionary& ofDict);
55 
57  template<class classType>
58  static label listDeleteVal(
59  List<classType>& listIn,
60  const classType valDel);
61 
63  static void writeMatrixBinary(
64  const Mat matIn,
65  const word prefix);
66 
68  static void writeMatrixASCII(
69  const Mat matIn,
70  const word prefix);
71 
73  static void readMatrixBinary(
74  Mat matIn,
75  const word prefix);
76 
78  static void writeVectorASCII(
79  const Vec vecIn,
80  const word prefix);
81 
83  static void readVectorBinary(
84  Vec vecIn,
85  const word prefix);
86 
88  static void writeVectorBinary(
89  const Vec vecIn,
90  const word prefix);
91 
93  static void boundVar(
94  const dictionary& allOptions,
95  volScalarField& var,
96  const label printToScreen);
97 
99  static void boundVar(
100  const dictionary& allOptions,
101  volVectorField& var,
102  const label printToScreen);
103 
105  static label isValueCloseToRef(
106  const scalar val,
107  const scalar refVal,
108  const scalar tol = 1.0e-6);
109 
111  static globalIndex genGlobalIndex(const label localIndexSize);
112 
115 
117  static void* pyCalcBeta;
119 
120  static void* pyCalcBetaJacVecProd;
122 
123  static void* pySetModelName;
125 
127  static void primalResidualControl(
128  const SolverPerformance<scalar>& solverP,
129  const label printToScreen,
130  const word varName,
131  scalar& primalMaxRes);
132 
134  static void primalResidualControl(
135  const SolverPerformance<vector>& solverP,
136  const label printToScreen,
137  const word varName,
138  scalar& primalMaxRes);
139 
140  static label myFindCell(
141  const primitiveMesh& mesh,
142  const point& point);
143 
144  static label isFieldReadable(
145  const fvMesh& mesh,
146  const word fieldName,
147  const word fieldType);
148 
150  template<class classType>
151  static void swapLists(List<classType>& a, List<classType>& b);
152 };
153 
154 template<class classType>
156  List<classType>& listIn,
157  const classType valDel)
158 {
159  /*
160  Description:
161  Delete a value in the list, the sequence of other
162  elements will be preserved
163 
164  Input:
165  listIn: the list to delete value
166  valDel: value to delete
167 
168  Output:
169  foundVal: 1 means the requested value is found and deleted
170 
171  Example:
172  List<word> listIn={"apple","orange"};
173  listDeleteVal<word>(listIn,"apple");
174  Now listIn will be {orange"}
175  NOTE: if list has multiple val, they will be all deleted
176  */
177 
178  // we will create a new list and delete the value
179  List<classType> listNew;
180 
181  label foundVal = 0;
182  forAll(listIn, idxI)
183  {
184  const classType& val = listIn[idxI];
185  if (val == valDel)
186  {
187  foundVal = 1;
188  // do not append
189  }
190  else
191  {
192  listNew.append(val);
193  }
194  }
195 
196  listIn.clear();
197  listIn = listNew;
198 
199  return foundVal;
200 }
201 
203 template<class classType>
204 void DAUtility::swapLists(List<classType>& a, List<classType>& b)
205 {
206  List<classType> temp = a;
207  a = b;
208  b = temp;
209 }
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #endif
218 
219 // ************************************************************************* //
Foam::DAUtility::readMatrixBinary
static void readMatrixBinary(Mat matIn, const word prefix)
read petsc matrix in binary format
Definition: DAUtility.C:378
Foam::DAUtility::pyCalcBetaJacVecProdInterface
static pyJacVecProdInterface pyCalcBetaJacVecProdInterface
Definition: DAUtility.H:121
Foam::DAUtility::readVectorBinary
static void readVectorBinary(Vec vecIn, const word prefix)
read petsc vector in binary format
Definition: DAUtility.C:282
Foam::DAUtility::pyCalcBeta
static void * pyCalcBeta
define a function pointer template for Python call back
Definition: DAUtility.H:117
Foam::DAUtility::DAUtility
DAUtility()
Constructors.
Definition: DAUtility.C:16
Foam::DAUtility::pyCalcBetaInterface
static pyComputeInterface pyCalcBetaInterface
Definition: DAUtility.H:118
Foam::DAUtility::~DAUtility
virtual ~DAUtility()
Destructor.
Definition: DAUtility.C:20
IOMRFZoneListDF.H
Foam::pySetCharInterface
void(* pySetCharInterface)(const char *, void *)
Definition: DAUtility.H:29
Foam::DAUtility
Definition: DAUtility.H:34
Foam::DAUtility::pyCalcBetaJacVecProd
static void * pyCalcBetaJacVecProd
Definition: DAUtility.H:120
Foam::DAUtility::pySetModelNameInterface
static pySetCharInterface pySetModelNameInterface
Definition: DAUtility.H:124
Foam::DAUtility::angleOfAttackRadForwardAD
static scalar angleOfAttackRadForwardAD
angle of attack in radian used in forward mode AD
Definition: DAUtility.H:114
Foam::DAUtility::isValueCloseToRef
static label isValueCloseToRef(const scalar val, const scalar refVal, const scalar tol=1.0e-6)
check whether a value is close to a reference value by a tolerance
Definition: DAUtility.C:707
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAUtility::pyDict2OFDict
static void pyDict2OFDict(PyObject *pyDict, dictionary &ofDict)
convert a python dictionary object to OpenFoam dictionary
Definition: DAUtility.C:24
Foam::DAUtility::pySetModelName
static void * pySetModelName
Definition: DAUtility.H:123
Foam::DAUtility::boundVar
static void boundVar(const dictionary &allOptions, volScalarField &var, const label printToScreen)
bound a volScalar variable based on parametes defined in DAOption::allOptions_
Definition: DAUtility.C:475
Foam::DAUtility::writeVectorBinary
static void writeVectorBinary(const Vec vecIn, const word prefix)
write petsc vector in binary format
Definition: DAUtility.C:315
Foam::DAUtility::swapLists
static void swapLists(List< classType > &a, List< classType > &b)
swap two lists
Definition: DAUtility.H:204
Foam::DAUtility::primalResidualControl
static void primalResidualControl(const SolverPerformance< scalar > &solverP, const label printToScreen, const word varName, scalar &primalMaxRes)
control when to print the residual and also compute the maxInitRes
Definition: DAUtility.C:734
Foam
Definition: checkGeometry.C:32
Foam::pyJacVecProdInterface
void(* pyJacVecProdInterface)(const double *, double *, int, const double *, const double *, int, void *)
Definition: DAUtility.H:28
allOptions
const dictionary & allOptions
Definition: createRefsPimpleDyM.H:14
Foam::DAUtility::listDeleteVal
static label listDeleteVal(List< classType > &listIn, const classType valDel)
delete a value in the list
Definition: DAUtility.H:155
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAUtility::writeMatrixBinary
static void writeMatrixBinary(const Mat matIn, const word prefix)
write petsc matrix in binary format
Definition: DAUtility.C:411
Foam::DAUtility::myFindCell
static label myFindCell(const primitiveMesh &mesh, const point &point)
Definition: DAUtility.C:803
Foam::DAUtility::writeVectorASCII
static void writeVectorASCII(const Vec vecIn, const word prefix)
write petsc vector in ascii format
Definition: DAUtility.C:347
Foam::DAUtility::writeMatrixASCII
static void writeMatrixASCII(const Mat matIn, const word prefix)
write petsc matrix in ascii format
Definition: DAUtility.C:443
Foam::pyComputeInterface
void(* pyComputeInterface)(const double *, int, double *, int, void *)
Definition: DAUtility.H:27
Foam::DAUtility::isFieldReadable
static label isFieldReadable(const fvMesh &mesh, const word fieldName, const word fieldType)
Definition: DAUtility.C:817
Foam::DAUtility::genGlobalIndex
static globalIndex genGlobalIndex(const label localIndexSize)
generate global index numbering for local-global index transferring
Definition: DAUtility.C:655