pyDASolvers.pyx
Go to the documentation of this file.
1 
2 # distutils: language = c++
3 # distutils: sources = DASolvers.C
4 
5 """
6  DAFoam : Discrete Adjoint with OpenFOAM
7  Version : v3
8 
9  Description:
10  Cython wrapper functions that call OpenFOAM libraries defined
11  in the *.C and *.H files. The python naming convention is to
12  add "py" before the C++ class name
13 """
14 
15 # for using Petsc
16 from petsc4py.PETSc cimport Vec, PetscVec, Mat, PetscMat, KSP, PetscKSP
17 cimport numpy as np
18 np.import_array() # initialize C API to call PyArray_SimpleNewFromData
19 
20 cdef public api CPointerToPyArray(const double* data, int size) with gil:
21  if not (data and size >= 0): raise ValueError
22  cdef np.npy_intp dims = size
23  return np.PyArray_SimpleNewFromData(1, &dims, np.NPY_DOUBLE, <void*>data)
24 
25 ctypedef void (*pyComputeInterface)(const double *, int, double *, int, void *)
26 ctypedef void (*pyJacVecProdInterface)(const double *, double *, int, const double *, const double *, int, void *)
27 
28 cdef void pyCalcBetaCallBack(const double* inputs, int n, double* outputs, int m, void *func):
29  inputs_data = CPointerToPyArray(inputs, n)
30  outputs_data = CPointerToPyArray(outputs, m)
31  (<object>func)(inputs_data, n, outputs_data, m)
32 
33 cdef void pyCalcBetaJacVecProdCallBack(const double* inputs, double* inputs_b, int n, const double* outputs, const double* outputs_b, int m, void *func):
34  inputs_data = CPointerToPyArray(inputs, n)
35  inputs_b_data = CPointerToPyArray(inputs_b, n)
36  outputs_data = CPointerToPyArray(outputs, m)
37  outputs_b_data = CPointerToPyArray(outputs_b, m)
38  (<object>func)(inputs_data, inputs_b_data, n, outputs_data, outputs_b_data, m)
39 
40 # declare cpp functions
41 cdef extern from "DASolvers.H" namespace "Foam":
42  cppclass DASolvers:
43  DASolvers(char *, object) except +
44  void initSolver()
45  int solvePrimal(PetscVec, PetscVec)
46  void calcdRdWT(PetscVec, PetscVec, int, PetscMat)
47  void calcdRdWTPsiAD(PetscVec, PetscVec, PetscVec, PetscVec)
48  void initializedRdWTMatrixFree(PetscVec, PetscVec)
49  void destroydRdWTMatrixFree()
50  void calcdFdW(PetscVec, PetscVec, char *, PetscVec)
51  void calcdFdWAD(PetscVec, PetscVec, char *, PetscVec)
52  void createMLRKSP(PetscMat, PetscMat, PetscKSP)
53  void createMLRKSPMatrixFree(PetscMat, PetscKSP)
54  void solveLinearEqn(PetscKSP, PetscVec, PetscVec)
55  void calcdRdBC(PetscVec, PetscVec, char *, PetscMat)
56  void calcdFdBC(PetscVec, PetscVec, char *, char *, PetscVec)
57  void calcdFdBCAD(PetscVec, PetscVec, char *, char *, PetscVec)
58  void calcdRdAOA(PetscVec, PetscVec, char *, PetscMat)
59  void calcdFdAOA(PetscVec, PetscVec, char *, char *, PetscVec)
60  void calcdRdFFD(PetscVec, PetscVec, char *, PetscMat)
61  void calcdRdXvTPsiAD(PetscVec, PetscVec, PetscVec, PetscVec)
62  void calcdForcedXvAD(PetscVec, PetscVec, PetscVec, PetscVec)
63  void calcdAcousticsdXvAD(PetscVec, PetscVec, PetscVec, PetscVec, char*, char*)
64  void calcdRdActTPsiAD(PetscVec, PetscVec, PetscVec, char*, PetscVec)
65  void calcdForcedWAD(PetscVec, PetscVec, PetscVec, PetscVec)
66  void calcdAcousticsdWAD(PetscVec, PetscVec, PetscVec, PetscVec, char*, char*)
67  void calcdFdACT(PetscVec, PetscVec, char *, char*, char*, PetscVec)
68  void calcdFdACTAD(PetscVec, PetscVec, char *, char*, PetscVec)
69  void calcdRdAOATPsiAD(PetscVec, PetscVec, PetscVec, char*, PetscVec)
70  void calcdRdBCTPsiAD(PetscVec, PetscVec, PetscVec, char*, PetscVec)
71  void calcdFdFFD(PetscVec, PetscVec, char *, char *, PetscVec)
72  void calcdFdXvAD(PetscVec, PetscVec, char *, char*, PetscVec)
73  void calcdRdACT(PetscVec, PetscVec, char *, char *, PetscMat)
74  void calcdRdFieldTPsiAD(PetscVec, PetscVec, PetscVec, char *, PetscVec)
75  void calcdFdFieldAD(PetscVec, PetscVec, char *, char *, PetscVec)
76  void calcdRdThermalTPsiAD(double *, double *, double *, double *, double *)
77  void calcdRdWOldTPsiAD(int, PetscVec, PetscVec)
78  void convertMPIVec2SeqVec(PetscVec, PetscVec)
79  void syncDAOptionToActuatorDVs()
80  void updateOFField(PetscVec)
81  void updateOFMesh(PetscVec)
82  void setdXvdFFDMat(PetscMat)
83  void setFFD2XvSeedVec(PetscVec)
84  int getGlobalXvIndex(int, int)
85  void ofField2StateVec(PetscVec)
86  void stateVec2OFField(PetscVec)
87  int getNLocalAdjointStates()
88  int getNLocalAdjointBoundaryStates()
89  int getNLocalCells()
90  int getNLocalPoints()
91  int getNCouplingFaces()
92  int getNCouplingPoints()
93  int checkMesh()
94  double getObjFuncValue(char *)
95  void calcCouplingFaceCoords(double *, double *)
96  void calcCouplingFaceCoordsAD(double *, double *, double *)
97  void getForces(PetscVec, PetscVec, PetscVec)
98  void getThermal(double *, double *, double *)
99  void getThermalAD(char *, double *, double *, double *, double *)
100  void setThermal(double *)
101  void getAcousticData(PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, char*)
102  void printAllOptions()
103  void updateDAOption(object)
104  double getPrevPrimalSolTime()
105  # functions for unit tests
106  void pointVec2OFMesh(PetscVec)
107  void ofMesh2PointVec(PetscVec)
108  void resVec2OFResField(PetscVec)
109  void ofResField2ResVec(PetscVec)
110  void writeMatrixBinary(PetscMat, char *)
111  void writeMatrixASCII(PetscMat, char *)
112  void readMatrixBinary(PetscMat, char *)
113  void writeVectorASCII(PetscVec, char *)
114  void readVectorBinary(PetscVec, char *)
115  void writeVectorBinary(PetscVec, char *)
116  void setTimeInstanceField(int)
117  void setTimeInstanceVar(char *, PetscMat, PetscMat, PetscVec, PetscVec)
118  double getTimeInstanceObjFunc(int, char *)
119  void setFieldValue4GlobalCellI(char *, double, int, int)
120  void setFieldValue4LocalCellI(char *, double, int, int)
121  void updateBoundaryConditions(char *, char *)
122  void calcPrimalResidualStatistics(char *)
123  double getForwardADDerivVal(char *)
124  void calcResidualVec(PetscVec)
125  void setPrimalBoundaryConditions(int)
126  void calcFvSource(char *, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec)
127  void calcdFvSourcedInputsTPsiAD(char *, char *, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec, PetscVec)
128  void calcForceProfile(PetscVec, PetscVec, PetscVec, PetscVec)
129  void calcdForcedStateTPsiAD(char *, PetscVec, PetscVec, PetscVec, PetscVec)
130  int runFPAdj(PetscVec, PetscVec, PetscVec, PetscVec)
131  void initTensorFlowFuncs(pyComputeInterface, void *, pyJacVecProdInterface, void *)
132 
133 # create python wrappers that call cpp functions
134 cdef class pyDASolvers:
135 
136  # define a class pointer for cpp functions
137  cdef:
138  DASolvers * _thisptr
139 
140  # initialize this class pointer with NULL
141  def __cinit__(self):
142  self._thisptr = NULL
143 
144  # deallocate the class pointer, and
145  # make sure we don't have memory leak
146  def __dealloc__(self):
147  if self._thisptr != NULL:
148  del self._thisptr
149 
150  # point the class pointer to the cpp class constructor
151  def __init__(self, argsAll, pyOptions):
152  """
153  Parameters
154  ----------
155 
156  argsAll: char
157  Chars that contains all the arguments
158  for running OpenFOAM solvers, including
159  the name of the solver.
160 
161  pyOptions: dict
162  Dictionary that defines all the options
163  in DAFoam
164 
165  Examples
166  --------
167  solver = pyDASolvers("DASolvers -parallel -python", aeroOptions)
168  """
169  self._thisptr = new DASolvers(argsAll, pyOptions)
170 
171  # wrap all the other member functions in the cpp class
172  def initSolver(self):
173  self._thisptr.initSolver()
174 
175  def solvePrimal(self, Vec xvVec, Vec wVec):
176  return self._thisptr.solvePrimal(xvVec.vec, wVec.vec)
177 
178  def calcdRdWT(self, Vec xvVec, Vec wVec, isPC, Mat dRdWT):
179  self._thisptr.calcdRdWT(xvVec.vec, wVec.vec, isPC, dRdWT.mat)
180 
181  def calcdRdWTPsiAD(self, Vec xvVec, Vec wVec, Vec psi, Vec dRdWTPsi):
182  self._thisptr.calcdRdWTPsiAD(xvVec.vec, wVec.vec, psi.vec, dRdWTPsi.vec)
183 
184  def initializedRdWTMatrixFree(self, Vec xvVec, Vec wVec):
185  self._thisptr.initializedRdWTMatrixFree(xvVec.vec, wVec.vec)
186 
187  def destroydRdWTMatrixFree(self):
188  self._thisptr.destroydRdWTMatrixFree()
189 
190  def calcdFdW(self, Vec xvVec, Vec wVec, objFuncName, Vec dFdW):
191  self._thisptr.calcdFdW(xvVec.vec, wVec.vec, objFuncName, dFdW.vec)
192 
193  def calcdFdWAD(self, Vec xvVec, Vec wVec, objFuncName, Vec dFdW):
194  self._thisptr.calcdFdWAD(xvVec.vec, wVec.vec, objFuncName, dFdW.vec)
195 
196  def createMLRKSP(self, Mat jacMat, Mat jacPCMat, KSP myKSP):
197  self._thisptr.createMLRKSP(jacMat.mat, jacPCMat.mat, myKSP.ksp)
198 
199  def createMLRKSPMatrixFree(self, Mat jacPCMat, KSP myKSP):
200  self._thisptr.createMLRKSPMatrixFree(jacPCMat.mat, myKSP.ksp)
201 
202  def solveLinearEqn(self, KSP myKSP, Vec rhsVec, Vec solVec):
203  self._thisptr.solveLinearEqn(myKSP.ksp, rhsVec.vec, solVec.vec)
204 
205  def calcdRdBC(self, Vec xvVec, Vec wVec, designVarName, Mat dRdBC):
206  self._thisptr.calcdRdBC(xvVec.vec, wVec.vec, designVarName, dRdBC.mat)
207 
208  def calcdFdBC(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdBC):
209  self._thisptr.calcdFdBC(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdBC.vec)
210 
211  def calcdFdBCAD(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdBC):
212  self._thisptr.calcdFdBCAD(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdBC.vec)
213 
214  def calcdRdAOA(self, Vec xvVec, Vec wVec, designVarName, Mat dRdAOA):
215  self._thisptr.calcdRdAOA(xvVec.vec, wVec.vec, designVarName, dRdAOA.mat)
216 
217  def calcdFdAOA(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdAOA):
218  self._thisptr.calcdFdAOA(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdAOA.vec)
219 
220  def calcdRdFFD(self, Vec xvVec, Vec wVec, designVarName, Mat dRdFFD):
221  self._thisptr.calcdRdFFD(xvVec.vec, wVec.vec, designVarName, dRdFFD.mat)
222 
223  def calcdRdXvTPsiAD(self, Vec xvVec, Vec wVec, Vec psi, Vec dRdXvTPsi):
224  self._thisptr.calcdRdXvTPsiAD(xvVec.vec, wVec.vec, psi.vec, dRdXvTPsi.vec)
225 
226  def calcdForcedXvAD(self, Vec xvVec, Vec wVec, Vec fBarVec, Vec dForcedXv):
227  self._thisptr.calcdForcedXvAD(xvVec.vec, wVec.vec, fBarVec.vec, dForcedXv.vec)
228 
229  def calcdAcousticsdXvAD(self, Vec xvVec, Vec wVec, Vec fBarVec, Vec dAcoudXv, varName, groupName):
230  self._thisptr.calcdAcousticsdXvAD(xvVec.vec, wVec.vec, fBarVec.vec, dAcoudXv.vec, varName, groupName)
231 
232  def calcdRdActTPsiAD(self, Vec xvVec, Vec wVec, Vec psi, designVarName, Vec dRdActTPsi):
233  self._thisptr.calcdRdActTPsiAD(xvVec.vec, wVec.vec, psi.vec, designVarName, dRdActTPsi.vec)
234 
235  def calcdForcedWAD(self, Vec xvVec, Vec wVec, Vec fBarVec, Vec dForcedW):
236  self._thisptr.calcdForcedWAD(xvVec.vec, wVec.vec, fBarVec.vec, dForcedW.vec)
237 
238  def calcdAcousticsdWAD(self, Vec xvVec, Vec wVec, Vec fBarVec, Vec dAcoudW, varName, groupName):
239  self._thisptr.calcdAcousticsdWAD(xvVec.vec, wVec.vec, fBarVec.vec, dAcoudW.vec, varName, groupName)
240 
241  def calcdFdACTAD(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdACT):
242  self._thisptr.calcdFdACTAD(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdACT.vec)
243 
244  def calcdFdACT(self, Vec xvVec, Vec wVec, objFuncName, designVarName, designVarType, Vec dFdACT):
245  self._thisptr.calcdFdACT(xvVec.vec, wVec.vec, objFuncName, designVarName, designVarType, dFdACT.vec)
246 
247  def calcdRdAOATPsiAD(self, Vec xvVec, Vec wVec, Vec psi, designVarName, Vec dRdAOATPsi):
248  self._thisptr.calcdRdAOATPsiAD(xvVec.vec, wVec.vec, psi.vec, designVarName, dRdAOATPsi.vec)
249 
250  def calcdRdBCTPsiAD(self, Vec xvVec, Vec wVec, Vec psi, designVarName, Vec dRdBCTPsi):
251  self._thisptr.calcdRdBCTPsiAD(xvVec.vec, wVec.vec, psi.vec, designVarName, dRdBCTPsi.vec)
252 
253  def calcdFdFFD(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdFFD):
254  self._thisptr.calcdFdFFD(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdFFD.vec)
255 
256  def calcdFdXvAD(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdXv):
257  self._thisptr.calcdFdXvAD(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdXv.vec)
258 
259  def calcdRdACT(self, Vec xvVec, Vec wVec, designVarName, designVarType, Mat dRdACT):
260  self._thisptr.calcdRdACT(xvVec.vec, wVec.vec, designVarName, designVarType, dRdACT.mat)
261 
262  def calcdRdFieldTPsiAD(self, Vec xvVec, Vec wVec, Vec psiVec, designVarName, Vec dRdFieldTPsi):
263  self._thisptr.calcdRdFieldTPsiAD(xvVec.vec, wVec.vec, psiVec.vec, designVarName, dRdFieldTPsi.vec)
264 
265  def calcdFdFieldAD(self, Vec xvVec, Vec wVec, objFuncName, designVarName, Vec dFdField):
266  self._thisptr.calcdFdFieldAD(xvVec.vec, wVec.vec, objFuncName, designVarName, dFdField.vec)
267 
268  def calcdRdThermalTPsiAD(self,
269  np.ndarray[double, ndim=1, mode="c"] volCoords,
270  np.ndarray[double, ndim=1, mode="c"] states,
271  np.ndarray[double, ndim=1, mode="c"] thermal,
272  np.ndarray[double, ndim=1, mode="c"] seeds,
273  np.ndarray[double, ndim=1, mode="c"] product):
274 
275  assert len(volCoords) == self.getNLocalPoints() * 3, "invalid array size!"
276  assert len(states) == self.getNLocalAdjointStates(), "invalid array size!"
277  assert len(thermal) == self.getNCouplingFaces() * 2, "invalid array size!"
278  assert len(seeds) == self.getNLocalAdjointStates(), "invalid array size!"
279  assert len(product) == self.getNCouplingFaces() * 2, "invalid array size!"
280 
281  cdef double *volCoords_data = <double*>volCoords.data
282  cdef double *states_data = <double*>states.data
283  cdef double *thermal_data = <double*>thermal.data
284  cdef double *seeds_data = <double*>seeds.data
285  cdef double *product_data = <double*>product.data
286 
287  self._thisptr.calcdRdThermalTPsiAD(
288  volCoords_data,
289  states_data,
290  thermal_data,
291  seeds_data,
292  product_data)
293 
294  def calcdRdWOldTPsiAD(self, oldTimeLevel, Vec psi, Vec dRdWOldTPsi):
295  self._thisptr.calcdRdWOldTPsiAD(oldTimeLevel, psi.vec, dRdWOldTPsi.vec)
296 
297  def convertMPIVec2SeqVec(self, Vec mpiVec, Vec seqVec):
298  self._thisptr.convertMPIVec2SeqVec(mpiVec.vec, seqVec.vec)
299 
300  def syncDAOptionToActuatorDVs(self):
301  self._thisptr.syncDAOptionToActuatorDVs()
302 
303  def updateOFField(self, Vec wVec):
304  self._thisptr.updateOFField(wVec.vec)
305 
306  def updateOFMesh(self, Vec xvVec):
307  self._thisptr.updateOFMesh(xvVec.vec)
308 
309  def setdXvdFFDMat(self, Mat dXvdFFDMat):
310  self._thisptr.setdXvdFFDMat(dXvdFFDMat.mat)
311 
312  def setFFD2XvSeedVec(self, Vec FFD2XvSeedVec):
313  self._thisptr.setFFD2XvSeedVec(FFD2XvSeedVec.vec)
314 
315  def getGlobalXvIndex(self, pointI, coordI):
316  return self._thisptr.getGlobalXvIndex(pointI, coordI)
317 
318  def ofField2StateVec(self, Vec stateVec):
319  self._thisptr.ofField2StateVec(stateVec.vec)
320 
321  def stateVec2OFField(self, Vec stateVec):
322  self._thisptr.stateVec2OFField(stateVec.vec)
323 
324  def pointVec2OFMesh(self, Vec xvVec):
325  self._thisptr.pointVec2OFMesh(xvVec.vec)
326 
327  def ofMesh2PointVec(self, Vec xvVec):
328  self._thisptr.ofMesh2PointVec(xvVec.vec)
329 
330  def resVec2OFResField(self, Vec rVec):
331  self._thisptr.resVec2OFResField(rVec.vec)
332 
333  def ofResField2ResVec(self, Vec rVec):
334  self._thisptr.ofResField2ResVec(rVec.vec)
335 
336  def getNLocalAdjointStates(self):
337  return self._thisptr.getNLocalAdjointStates()
338 
339  def getNCouplingFaces(self):
340  return self._thisptr.getNCouplingFaces()
341 
342  def getNCouplingPoints(self):
343  return self._thisptr.getNCouplingPoints()
344 
345  def getNLocalAdjointBoundaryStates(self):
346  return self._thisptr.getNLocalAdjointBoundaryStates()
347 
348  def getNLocalCells(self):
349  return self._thisptr.getNLocalCells()
350 
351  def getNLocalPoints(self):
352  return self._thisptr.getNLocalPoints()
353 
354  def checkMesh(self):
355  return self._thisptr.checkMesh()
356 
357  def getObjFuncValue(self, objFuncName):
358  return self._thisptr.getObjFuncValue(objFuncName)
359 
360  def calcCouplingFaceCoords(self,
361  np.ndarray[double, ndim=1, mode="c"] volCoords,
362  np.ndarray[double, ndim=1, mode="c"] surfCoords):
363 
364  assert len(volCoords) == self.getNLocalPoints() * 3, "invalid array size!"
365  assert len(surfCoords) == self.getNCouplingFaces() * 6, "invalid array size!"
366 
367  cdef double *volCoords_data = <double*>volCoords.data
368  cdef double *surfCoords_data = <double*>surfCoords.data
369 
370  self._thisptr.calcCouplingFaceCoords(volCoords_data, surfCoords_data)
371 
372  def calcCouplingFaceCoordsAD(self,
373  np.ndarray[double, ndim=1, mode="c"] volCoords,
374  np.ndarray[double, ndim=1, mode="c"] seeds,
375  np.ndarray[double, ndim=1, mode="c"] product):
376 
377  assert len(volCoords) == self.getNLocalPoints() * 3, "invalid array size!"
378  assert len(seeds) == self.getNCouplingFaces() * 6, "invalid array size!"
379  assert len(product) == self.getNLocalPoints() * 3, "invalid array size!"
380 
381  cdef double *volCoords_data = <double*>volCoords.data
382  cdef double *seeds_data = <double*>seeds.data
383  cdef double *product_data = <double*>product.data
384 
385  self._thisptr.calcCouplingFaceCoordsAD(volCoords_data, seeds_data, product_data)
386 
387  def getForces(self, Vec fX, Vec fY, Vec fZ):
388  self._thisptr.getForces(fX.vec, fY.vec, fZ.vec)
389 
390  def getThermal(self,
391  np.ndarray[double, ndim=1, mode="c"] volCoords,
392  np.ndarray[double, ndim=1, mode="c"] states,
393  np.ndarray[double, ndim=1, mode="c"] thermal):
394 
395  assert len(volCoords) == self.getNLocalPoints() * 3, "invalid array size!"
396  assert len(states) == self.getNLocalAdjointStates(), "invalid array size!"
397  assert len(thermal) == self.getNCouplingFaces() * 2, "invalid array size!"
398 
399  cdef double *volCoords_data = <double*>volCoords.data
400  cdef double *states_data = <double*>states.data
401  cdef double *thermal_data = <double*>thermal.data
402 
403  self._thisptr.getThermal(volCoords_data, states_data, thermal_data)
404 
405  def getThermalAD(self,
406  inputName,
407  np.ndarray[double, ndim=1, mode="c"] volCoords,
408  np.ndarray[double, ndim=1, mode="c"] states,
409  np.ndarray[double, ndim=1, mode="c"] seeds,
410  np.ndarray[double, ndim=1, mode="c"] product):
411 
412  assert len(volCoords) == self.getNLocalPoints() * 3, "invalid array size!"
413  assert len(states) == self.getNLocalAdjointStates(), "invalid array size!"
414  assert len(seeds) == self.getNCouplingFaces() * 2, "invalid array size!"
415  if inputName == "volCoords":
416  assert len(product) == self.getNLocalPoints() * 3, "invalid array size!"
417  elif inputName == "states":
418  assert len(product) == self.getNLocalAdjointStates(), "invalid array size!"
419  else:
420  print("invalid inputName!")
421  exit(1)
422 
423  cdef double *volCoords_data = <double*>volCoords.data
424  cdef double *states_data = <double*>states.data
425  cdef double *seeds_data = <double*>seeds.data
426  cdef double *product_data = <double*>product.data
427 
428  self._thisptr.getThermalAD(
429  inputName.encode(),
430  volCoords_data,
431  states_data,
432  seeds_data,
433  product_data)
434 
435  def setThermal(self,
436  np.ndarray[double, ndim=1, mode="c"] thermal):
437 
438  assert len(thermal) == self.getNCouplingFaces() * 2, "invalid array size!"
439  cdef double *thermal_data = <double*>thermal.data
440  self._thisptr.setThermal(thermal_data)
441 
442  def getAcousticData(self, Vec x, Vec y, Vec z, Vec nX, Vec nY, Vec nZ, Vec a, Vec fX, Vec fY, Vec fZ, groupName):
443  self._thisptr.getAcousticData(x.vec, y.vec, z.vec, nX.vec, nY.vec, nZ.vec, a.vec, fX.vec, fY.vec, fZ.vec, groupName)
444 
445  def getAcousticData(self, Vec x, Vec y, Vec z, Vec nX, Vec nY, Vec nZ, Vec a, Vec fX, Vec fY, Vec fZ, groupName):
446  self._thisptr.getAcousticData(x.vec, y.vec, z.vec, nX.vec, nY.vec, nZ.vec, a.vec, fX.vec, fY.vec, fZ.vec, groupName)
447 
448  def printAllOptions(self):
449  self._thisptr.printAllOptions()
450 
451  def updateDAOption(self, pyOptions):
452  self._thisptr.updateDAOption(pyOptions)
453 
454  def getPrevPrimalSolTime(self):
455  return self._thisptr.getPrevPrimalSolTime()
456 
457  def writeMatrixBinary(self, Mat magIn, prefix):
458  self._thisptr.writeMatrixBinary(magIn.mat, prefix)
459 
460  def writeMatrixASCII(self, Mat magIn, prefix):
461  self._thisptr.writeMatrixASCII(magIn.mat, prefix)
462 
463  def readMatrixBinary(self, Mat magIn, prefix):
464  self._thisptr.readMatrixBinary(magIn.mat, prefix)
465 
466  def writeVectorASCII(self, Vec vecIn, prefix):
467  self._thisptr.writeVectorASCII(vecIn.vec, prefix)
468 
469  def readVectorBinary(self, Vec vecIn, prefix):
470  self._thisptr.readVectorBinary(vecIn.vec, prefix)
471 
472  def writeVectorBinary(self, Vec vecIn, prefix):
473  self._thisptr.writeVectorBinary(vecIn.vec, prefix)
474 
475  def setTimeInstanceField(self, instanceI):
476  self._thisptr.setTimeInstanceField(instanceI)
477 
478  def setTimeInstanceVar(self, mode, Mat stateMat, Mat stateBCMat, Vec timeVec, Vec timeIdxVec):
479  self._thisptr.setTimeInstanceVar(mode, stateMat.mat, stateBCMat.mat, timeVec.vec, timeIdxVec.vec)
480 
481  def getTimeInstanceObjFunc(self, instanceI, objFuncName):
482  return self._thisptr.getTimeInstanceObjFunc(instanceI, objFuncName)
483 
484  def setFieldValue4GlobalCellI(self, fieldName, val, globalCellI, compI = 0):
485  self._thisptr.setFieldValue4GlobalCellI(fieldName, val, globalCellI, compI)
486 
487  def setFieldValue4LocalCellI(self, fieldName, val, localCellI, compI = 0):
488  self._thisptr.setFieldValue4LocalCellI(fieldName, val, localCellI, compI)
489 
490  def updateBoundaryConditions(self, fieldName, fieldType):
491  self._thisptr.updateBoundaryConditions(fieldName, fieldType)
492 
493  def calcPrimalResidualStatistics(self, mode):
494  self._thisptr.calcPrimalResidualStatistics(mode)
495 
496  def getForwardADDerivVal(self, objFuncName):
497  return self._thisptr.getForwardADDerivVal(objFuncName)
498 
499  def calcResidualVec(self, Vec resVec):
500  self._thisptr.calcResidualVec(resVec.vec)
501 
502  def setPrimalBoundaryConditions(self, printInfo):
503  self._thisptr.setPrimalBoundaryConditions(printInfo)
504 
505  def calcFvSource(self, propName, Vec aForce, Vec tForce, Vec rDist, Vec targetForce, Vec center, Vec fvSource):
506  self._thisptr.calcFvSource(propName, aForce.vec, tForce.vec, rDist.vec, targetForce.vec, center.vec, fvSource.vec)
507 
508  def calcdFvSourcedInputsTPsiAD(self, propName, mode, Vec aForce, Vec tForce, Vec rDist, Vec targetForce, Vec center, Vec psi, Vec dFvSource):
509  self._thisptr.calcdFvSourcedInputsTPsiAD(propName, mode, aForce.vec, tForce.vec, rDist.vec, targetForce.vec, center.vec, psi.vec, dFvSource.vec)
510 
511  def calcForceProfile(self, Vec center, Vec aForceL, Vec tForceL, Vec rDistL):
512  self._thisptr.calcForceProfile(center.vec, aForceL.vec, tForceL.vec, rDistL.vec)
513 
514  def calcdForcedStateTPsiAD(self, mode, Vec xv, Vec state, Vec psi, Vec prod):
515  self._thisptr.calcdForcedStateTPsiAD(mode, xv.vec, state.vec, psi.vec, prod.vec)
516 
517  def runFPAdj(self, Vec xvVec, Vec wVec, Vec dFdW, Vec psi):
518  return self._thisptr.runFPAdj(xvVec.vec, wVec.vec, dFdW.vec, psi.vec)
519 
520  def initTensorFlowFuncs(self, compute, jacVecProd):
521  self._thisptr.initTensorFlowFuncs(pyCalcBetaCallBack, <void*>compute, pyCalcBetaJacVecProdCallBack, <void*>jacVecProd)