DASolvers.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  Python wrapper library for DASolver
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DASolvers_H
12 #define DASolvers_H
13 
14 #include <petscksp.h>
15 #include "Python.h"
16 #include "DASolver.H"
17 
18 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
19 
20 namespace Foam
21 {
22 
23 /*---------------------------------------------------------------------------*\
24  Class DASolvers Declaration
25 \*---------------------------------------------------------------------------*/
26 
27 class DASolvers
28 {
29 
30 private:
32  DASolvers(const DASolvers&);
33 
35  void operator=(const DASolvers&);
36 
38  char* argsAll_;
39 
41  PyObject* pyOptions_;
42 
44  autoPtr<DASolver> DASolverPtr_;
45 
46 public:
47  // Constructors
48 
50  DASolvers(
51  char* argsAll,
52  PyObject* pyOptions);
53 
55  virtual ~DASolvers();
56 
58  void initSolver()
59  {
60  DASolverPtr_->initSolver();
61  }
62 
64  label solvePrimal()
65  {
66  return DASolverPtr_->solvePrimal();
67  }
68 
69  label getInputSize(
70  const word inputName,
71  const word inputType)
72  {
73  return DASolverPtr_->getInputSize(inputName, inputType);
74  }
75 
77  const word outputName,
78  const word outputType)
79  {
80  return DASolverPtr_->getOutputSize(outputName, outputType);
81  }
82 
83  void calcOutput(
84  const word outputName,
85  const word outputType,
86  double* output)
87  {
88  DASolverPtr_->calcOutput(outputName, outputType, output);
89  }
90 
92  const word inputName,
93  const word inputType)
94  {
95  return DASolverPtr_->getInputDistributed(inputName, inputType);
96  }
97 
99  const word outputName,
100  const word outputType)
101  {
102  return DASolverPtr_->getOutputDistributed(outputName, outputType);
103  }
104 
107  const word inputName,
108  const word inputType,
109  const double* input,
110  const word outputName,
111  const word outputType,
112  const double* seed,
113  double* product)
114  {
115  DASolverPtr_->calcJacTVecProduct(
116  inputName,
117  inputType,
118  input,
119  outputName,
120  outputType,
121  seed,
122  product);
123  }
124 
126  const word inputName,
127  const word inputType,
128  const int inputSize,
129  const double* input,
130  const double* seed)
131  {
132  DASolverPtr_->setSolverInput(
133  inputName,
134  inputType,
135  inputSize,
136  input,
137  seed);
138  }
139 
141  void runColoring()
142  {
143  DASolverPtr_->runColoring();
144  }
145 
147  void calcdRdWT(
148  const label isPC,
149  Mat dRdWT)
150  {
151  DASolverPtr_->calcdRdWT(isPC, dRdWT);
152  }
153 
156  Mat PCMat,
157  KSP ksp)
158  {
159  DASolverPtr_->updateKSPPCMat(PCMat, ksp);
160  }
161 
164  const Mat jacPCMat,
165  KSP ksp)
166  {
167  DASolverPtr_->createMLRKSPMatrixFree(jacPCMat, ksp);
168  }
169 
172  {
173  DASolverPtr_->initializedRdWTMatrixFree();
174  }
175 
178  {
179  DASolverPtr_->destroydRdWTMatrixFree();
180  }
181 
184  const KSP ksp,
185  const Vec rhsVec,
186  Vec solVec)
187  {
188  return DASolverPtr_->solveLinearEqn(ksp, rhsVec, solVec);
189  }
190 
193  const label oldTimeLevel,
194  const double* psi,
195  double* dRdWOldTPsi)
196  {
197  DASolverPtr_->calcdRdWOldTPsiAD(oldTimeLevel, psi, dRdWOldTPsi);
198  }
199 
201  void updateOFFields(const double* states)
202  {
203 #if !defined(CODI_ADR) && !defined(CODI_ADF)
204  DASolverPtr_->updateOFFields(states);
205 #else
206  label localSize = this->getNLocalAdjointStates();
207  scalar* statesArray = new scalar[localSize];
208  for (label i = 0; i < localSize; i++)
209  {
210  statesArray[i] = states[i];
211  }
212  DASolverPtr_->updateOFFields(statesArray);
213  delete[] statesArray;
214 #endif
215  }
216 
218  void getOFMeshPoints(double* points)
219  {
220  DASolverPtr_->getOFMeshPoints(points);
221  }
222 
224  void getOFFields(double* states)
225  {
226 #if !defined(CODI_ADR) && !defined(CODI_ADF)
227  DASolverPtr_->getOFFields(states);
228 #else
229  label localSize = this->getNLocalAdjointStates();
230  scalar* statesArray = new scalar[localSize];
231  DASolverPtr_->getOFFields(statesArray);
232  for (label i = 0; i < localSize; i++)
233  {
234  states[i] = statesArray[i].value();
235  }
236  delete[] statesArray;
237 #endif
238  }
239 
242  const word fieldName,
243  const word fieldType,
244  double* field)
245  {
246  DASolverPtr_->getOFField(fieldName, fieldType, field);
247  }
248 
250  void updateOFMesh(const double* volCoords)
251  {
252 #if !defined(CODI_ADR) && !defined(CODI_ADF)
253  DASolverPtr_->updateOFMesh(volCoords);
254 #else
255  label localSize = this->getNLocalPoints() * 3;
256  scalar* volCoordsArray = new scalar[localSize];
257  for (label i = 0; i < localSize; i++)
258  {
259  volCoordsArray[i] = volCoords[i];
260  }
261  DASolverPtr_->updateOFMesh(volCoordsArray);
262  delete[] volCoordsArray;
263 #endif
264  }
265 
268  const label idxPoint,
269  const label idxCoord) const
270  {
271  return DASolverPtr_->getGlobalXvIndex(idxPoint, idxCoord);
272  }
273 
275  label checkMesh() const
276  {
277  return DASolverPtr_->checkMesh();
278  }
279 
280  // return the number of local mesh points
281  label getNLocalPoints() const
282  {
283  return DASolverPtr_->getNLocalPoints();
284  }
285 
288  {
289  return DASolverPtr_->getNLocalAdjointStates();
290  }
291 
294  {
295  return DASolverPtr_->getNLocalAdjointBoundaryStates();
296  }
297 
299  label getNLocalCells() const
300  {
301  return DASolverPtr_->getNLocalCells();
302  }
303 
305  double getTimeOpFuncVal(const word functionName)
306  {
307  return DASolverPtr_->getTimeOpFuncVal(functionName);
308  }
309 
311  double getdFScaling(const word functionName, const label timeIdx)
312  {
313  double returnVal = 0.0;
314  assignValueCheckAD(returnVal, DASolverPtr_->getdFScaling(functionName, timeIdx));
315  return returnVal;
316  }
317 
319  const double* volCoords,
320  double* surfCoords)
321  {
322 #if !defined(CODI_ADR) && !defined(CODI_ADF)
323  DASolverPtr_->calcCouplingFaceCoords(volCoords, surfCoords);
324  return;
325 #endif
326  FatalErrorIn("calcCouplingFaceCoords") << "should not call this func in AD mode!"
327  << abort(FatalError);
328  }
329 
332  {
333  double returnVal = 0.0;
334  assignValueCheckAD(returnVal, DASolverPtr_->getElapsedClockTime());
335  return returnVal;
336  }
337 
340  {
341  double returnVal = 0.0;
342  assignValueCheckAD(returnVal, DASolverPtr_->getElapsedCpuTime());
343  return returnVal;
344  }
345 
347  label getNRegressionParameters(word modelName)
348  {
349  return DASolverPtr_->getNRegressionParameters(modelName);
350  }
351 
354  {
355  DASolverPtr_->printAllOptions();
356  }
357 
360  {
361  return DASolverPtr_->hasVolCoordInput();
362  }
363 
365  void updateDAOption(PyObject* pyOptions)
366  {
367  DASolverPtr_->updateDAOption(pyOptions);
368  }
369 
372  {
373  double returnVal = 0.0;
374  assignValueCheckAD(returnVal, DASolverPtr_->getPrevPrimalSolTime());
375  return returnVal;
376  }
377 
379  {
380  DASolverPtr_->writeFailedMesh();
381  }
382 
385  scalar timeVal,
386  label timeLevel = 0)
387  {
388  DASolverPtr_->readStateVars(timeVal, timeLevel);
389  }
390 
392  void readMeshPoints(const scalar timeVal)
393  {
394  DASolverPtr_->readMeshPoints(timeVal);
395  }
396 
398  void writeMeshPoints(const double* points, const scalar timeVal)
399  {
400  DASolverPtr_->writeMeshPoints(points, timeVal);
401  }
402 
404  void calcPCMatWithFvMatrix(Mat PCMat, const label turbOnly)
405  {
406  DASolverPtr_->calcPCMatWithFvMatrix(PCMat, turbOnly);
407  }
408 
410  void setTime(scalar time, label timeIndex)
411  {
412  DASolverPtr_->setTime(time, timeIndex);
413  }
414 
417  {
418  return DASolverPtr_->getDdtSchemeOrder();
419  }
420 
422  double getEndTime()
423  {
424 #if !defined(CODI_ADR) && !defined(CODI_ADF)
425  return DASolverPtr_->getRunTime().endTime().value();
426 #endif
427  FatalErrorIn("getEndTime") << "should not call this func in AD mode!"
428  << abort(FatalError);
429  }
430 
432  double getDeltaT()
433  {
434 #if !defined(CODI_ADR) && !defined(CODI_ADF)
435  return DASolverPtr_->getRunTime().deltaT().value();
436 #endif
437  FatalErrorIn("getDeltaT") << "should not call this func in AD mode!"
438  << abort(FatalError);
439  }
440 
443  const word fieldName,
444  const word fieldType)
445  {
446  DASolverPtr_->updateBoundaryConditions(fieldName, fieldType);
447  }
448 
451  {
452  DASolverPtr_->updateStateBoundaryConditions();
453  }
454 
457  {
458  DASolverPtr_->calcPrimalResidualStatistics(mode, 0);
459  }
460 
461  void setPrimalBoundaryConditions(const label printInfo = 1)
462  {
463  DASolverPtr_->setPrimalBoundaryConditions(printInfo);
464  }
465 
466  label runFPAdj(
467  Vec dFdW,
468  Vec psi)
469  {
470  return DASolverPtr_->runFPAdj(dFdW, psi);
471  }
472 
474  Vec dFdW,
475  Vec psi)
476  {
477  return DASolverPtr_->solveAdjointFP(dFdW, psi);
478  }
479 
482  pyComputeInterface computeInterface,
483  void* compute,
484  pyJacVecProdInterface jacVecProdInterface,
485  void* jacVecProd,
486  pySetCharInterface setModelNameInterface,
487  void* setModelName)
488  {
489  DASolverPtr_->initTensorFlowFuncs(
490  computeInterface, compute, jacVecProdInterface, jacVecProd, setModelNameInterface, setModelName);
491  }
492 
495  const word function,
496  const double writeTime,
497  const double* psi)
498  {
499  DASolverPtr_->writeAdjointFields(function, writeTime, psi);
500  }
501 
504  const word name,
505  const double* dFdXs,
506  const double* Xs,
507  const label size,
508  const double timeName)
509  {
510  DASolverPtr_->writeSensMapSurface(name, dFdXs, Xs, size, timeName);
511  }
512 
515  const word name,
516  const double* dFdField,
517  const word fieldType,
518  const double timeName)
519  {
520  DASolverPtr_->writeSensMapField(name, dFdField, fieldType, timeName);
521  }
522 
524  double getLatestTime()
525  {
526  double latestTime = 0.0;
527 #if !defined(CODI_ADR) && !defined(CODI_ADF)
528  latestTime = DASolverPtr_->getLatestTime();
529 #else
530  scalar tmp = DASolverPtr_->getLatestTime();
531  latestTime = tmp.getValue();
532 #endif
533  return latestTime;
534  };
535 
538  {
539  DASolverPtr_->meanStatesToStates();
540  }
541 
544  {
545  DASolverPtr_->updateInputFieldUnsteady();
546  }
547 };
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 } // End namespace Foam
552 
553 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
554 
555 #endif
Foam::DASolvers
Definition: DASolvers.H:27
Foam::DASolvers::writeMeshPoints
void writeMeshPoints(const double *points, const scalar timeVal)
write the mesh points to the disk for the given timeVal
Definition: DASolvers.H:398
Foam::DASolvers::calcPrimalResidualStatistics
void calcPrimalResidualStatistics(const word mode)
Calculate the mean, max, and norm2 for all residuals and print it to screen.
Definition: DASolvers.H:456
Foam::DASolvers::getEndTime
double getEndTime()
return the endTime
Definition: DASolvers.H:422
Foam::DASolvers::updateBoundaryConditions
void updateBoundaryConditions(const word fieldName, const word fieldType)
update the boundary condition for a field
Definition: DASolvers.H:442
Foam::DASolvers::getGlobalXvIndex
label getGlobalXvIndex(const label idxPoint, const label idxCoord) const
basically, we call DAIndex::getGlobalXvIndex
Definition: DASolvers.H:267
Foam::DASolvers::getOFFields
void getOFFields(double *states)
Assign the OpenFOAM field values to the states array.
Definition: DASolvers.H:224
Foam::DASolvers::checkMesh
label checkMesh() const
basically, we call DASolver::checkMesh
Definition: DASolvers.H:275
Foam::DASolvers::getPrevPrimalSolTime
double getPrevPrimalSolTime()
get the solution time folder for previous primal solution
Definition: DASolvers.H:371
Foam::DASolvers::~DASolvers
virtual ~DASolvers()
Destructor.
Definition: DASolvers.C:25
Foam::DASolvers::writeAdjointFields
void writeAdjointFields(const word function, const double writeTime, const double *psi)
write the adjoint variables for all states
Definition: DASolvers.H:494
Foam::DASolvers::calcJacTVecProduct
void calcJacTVecProduct(const word inputName, const word inputType, const double *input, const word outputName, const word outputType, const double *seed, double *product)
calculate the Jacobian-matrix-transposed and vector product for product = [dOutput/dInput]^T * seed
Definition: DASolvers.H:106
Foam::DASolvers::createMLRKSPMatrixFree
void createMLRKSPMatrixFree(const Mat jacPCMat, KSP ksp)
create a multi-level, Richardson KSP object with matrix-free Jacobians
Definition: DASolvers.H:163
Foam::DASolvers::initializedRdWTMatrixFree
void initializedRdWTMatrixFree()
initialize matrix free dRdWT
Definition: DASolvers.H:171
Foam::DASolvers::getNRegressionParameters
label getNRegressionParameters(word modelName)
get the number of regression model parameters
Definition: DASolvers.H:347
Foam::DASolvers::getElapsedClockTime
double getElapsedClockTime()
return the elapsed clock time for testing speed
Definition: DASolvers.H:331
Foam::DASolvers::setSolverInput
void setSolverInput(const word inputName, const word inputType, const int inputSize, const double *input, const double *seed)
Definition: DASolvers.H:125
Foam::DASolvers::solveLinearEqn
label solveLinearEqn(const KSP ksp, const Vec rhsVec, Vec solVec)
solve the linear equation
Definition: DASolvers.H:183
Foam::pySetCharInterface
void(* pySetCharInterface)(const char *, void *)
Definition: DAUtility.H:29
Foam::DASolvers::getOFField
void getOFField(const word fieldName, const word fieldType, double *field)
get a field variable from OF layer
Definition: DASolvers.H:241
Foam::DASolvers::getOFMeshPoints
void getOFMeshPoints(double *points)
get the flatten mesh points coordinates
Definition: DASolvers.H:218
Foam::DASolvers::meanStatesToStates
void meanStatesToStates()
assign the mean states values to states
Definition: DASolvers.H:537
Foam::DASolvers::getElapsedCpuTime
double getElapsedCpuTime()
return the elapsed CPU time for testing speed
Definition: DASolvers.H:339
Foam::DASolvers::getDdtSchemeOrder
label getDdtSchemeOrder()
get the ddtScheme order
Definition: DASolvers.H:416
Foam::DASolvers::runColoring
void runColoring()
run the coloring solver
Definition: DASolvers.H:141
Foam::DASolvers::initTensorFlowFuncs
void initTensorFlowFuncs(pyComputeInterface computeInterface, void *compute, pyJacVecProdInterface jacVecProdInterface, void *jacVecProd, pySetCharInterface setModelNameInterface, void *setModelName)
initialize the call back function pointer
Definition: DASolvers.H:481
Foam::DASolvers::getOutputDistributed
label getOutputDistributed(const word outputName, const word outputType)
Definition: DASolvers.H:98
Foam::DASolvers::runFPAdj
label runFPAdj(Vec dFdW, Vec psi)
Definition: DASolvers.H:466
Foam::DASolvers::printAllOptions
void printAllOptions()
call DASolver::printAllOptions
Definition: DASolvers.H:353
Foam::DASolvers::readMeshPoints
void readMeshPoints(const scalar timeVal)
read the mesh points from the disk and run movePoints to deform the mesh
Definition: DASolvers.H:392
Foam::DASolvers::solveAdjointFP
label solveAdjointFP(Vec dFdW, Vec psi)
Definition: DASolvers.H:473
Foam::DASolvers::getLatestTime
double getLatestTime()
get the latest solution time
Definition: DASolvers.H:524
Foam::DASolvers::getNLocalPoints
label getNLocalPoints() const
Definition: DASolvers.H:281
Foam::DASolvers::writeFailedMesh
void writeFailedMesh()
Definition: DASolvers.H:378
Foam::DASolvers::updateInputFieldUnsteady
void updateInputFieldUnsteady()
assign the inputFieldUnsteady values to the OF field vars
Definition: DASolvers.H:543
Foam::DASolvers::calcCouplingFaceCoords
void calcCouplingFaceCoords(const double *volCoords, double *surfCoords)
Definition: DASolvers.H:318
Foam::DASolvers::getOutputSize
label getOutputSize(const word outputName, const word outputType)
Definition: DASolvers.H:76
Foam::DASolvers::updateOFFields
void updateOFFields(const double *states)
Update the OpenFOAM field values (including both internal and boundary fields) based on the states ar...
Definition: DASolvers.H:201
Foam::DASolvers::calcdRdWT
void calcdRdWT(const label isPC, Mat dRdWT)
compute dRdWT
Definition: DASolvers.H:147
Foam::DASolvers::getNLocalAdjointStates
label getNLocalAdjointStates() const
return the number of local adjoint states
Definition: DASolvers.H:287
Foam::DASolvers::writeSensMapField
void writeSensMapField(const word name, const double *dFdField, const word fieldType, const double timeName)
write the sensitivity map for the entire field
Definition: DASolvers.H:514
Foam::DASolvers::calcOutput
void calcOutput(const word outputName, const word outputType, double *output)
Definition: DASolvers.H:83
assignValueCheckAD
#define assignValueCheckAD(valA, valB)
Definition: DAMacroFunctions.H:103
Foam::DASolvers::getInputSize
label getInputSize(const word inputName, const word inputType)
Definition: DASolvers.H:69
Foam::DASolvers::updateOFMesh
void updateOFMesh(const double *volCoords)
Update the OpenFoam mesh point coordinates based on the point vector xvVec.
Definition: DASolvers.H:250
Foam::DASolvers::calcPCMatWithFvMatrix
void calcPCMatWithFvMatrix(Mat PCMat, const label turbOnly)
calculate the PC mat using fvMatrix
Definition: DASolvers.H:404
Foam::DASolvers::updateKSPPCMat
void updateKSPPCMat(Mat PCMat, KSP ksp)
Update the preconditioner matrix for the ksp object.
Definition: DASolvers.H:155
Foam::DASolvers::getdFScaling
double getdFScaling(const word functionName, const label timeIdx)
get the scaling factor for dF/d? computation
Definition: DASolvers.H:311
Foam::DASolvers::setPrimalBoundaryConditions
void setPrimalBoundaryConditions(const label printInfo=1)
Definition: DASolvers.H:461
Foam::DASolvers::getDeltaT
double getDeltaT()
return the deltaT
Definition: DASolvers.H:432
Foam
Definition: checkGeometry.C:32
Foam::DASolvers::updateStateBoundaryConditions
void updateStateBoundaryConditions()
update the boundary conditions for all states and intermediate variables
Definition: DASolvers.H:450
Foam::pyJacVecProdInterface
void(* pyJacVecProdInterface)(const double *, double *, int, const double *, const double *, int, void *)
Definition: DAUtility.H:28
Foam::DASolvers::destroydRdWTMatrixFree
void destroydRdWTMatrixFree()
destroy matrix free dRdWT
Definition: DASolvers.H:177
Foam::DASolvers::readStateVars
void readStateVars(scalar timeVal, label timeLevel=0)
read the state variables from the disk and assign the value to the prescribe time level
Definition: DASolvers.H:384
Foam::DASolvers::updateDAOption
void updateDAOption(PyObject *pyOptions)
update the allOptions_ dict in DAOption based on the pyOptions from pyDAFoam
Definition: DASolvers.H:365
Foam::DASolvers::initSolver
void initSolver()
initialize fields and variables
Definition: DASolvers.H:58
Foam::DASolvers::getNLocalAdjointBoundaryStates
label getNLocalAdjointBoundaryStates() const
return the number of local adjoint boundary states
Definition: DASolvers.H:293
Foam::DASolvers::setTime
void setTime(scalar time, label timeIndex)
setTime for OF fields
Definition: DASolvers.H:410
Foam::DASolvers::hasVolCoordInput
label hasVolCoordInput()
whether the volCoord input is defined
Definition: DASolvers.H:359
Foam::pyComputeInterface
void(* pyComputeInterface)(const double *, int, double *, int, void *)
Definition: DAUtility.H:27
Foam::DASolvers::solvePrimal
label solvePrimal()
solve the primal equations
Definition: DASolvers.H:64
DASolver.H
Foam::DASolvers::getNLocalCells
label getNLocalCells() const
return the number of local cells
Definition: DASolvers.H:299
Foam::DASolvers::writeSensMapSurface
void writeSensMapSurface(const word name, const double *dFdXs, const double *Xs, const label size, const double timeName)
write the sensitivity map for all wall surfaces
Definition: DASolvers.H:503
Foam::DASolvers::getTimeOpFuncVal
double getTimeOpFuncVal(const word functionName)
return the value of objective based on timeOp
Definition: DASolvers.H:305
Foam::DASolvers::getInputDistributed
label getInputDistributed(const word inputName, const word inputType)
Definition: DASolvers.H:91
Foam::DASolvers::calcdRdWOldTPsiAD
void calcdRdWOldTPsiAD(const label oldTimeLevel, const double *psi, double *dRdWOldTPsi)
compute dRdWOld^T*Psi
Definition: DASolvers.H:192
psi
const volScalarField & psi
Definition: createRefsRhoPimple.H:14
dafoam_plot3dtransform.mode
mode
Definition: dafoam_plot3dtransform.py:21