DAUtility.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAUtility.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 // Constructors
17 {
18 }
19 
21 {
22 }
23 
25  PyObject* pyDict,
26  dictionary& ofDict)
27 {
28  /*
29  Description:
30  Parse a Python dictionary to an OpenFOAM dictionary
31  Only support a certain number of data types
32 
33  Input:
34  pyDict: Pytion dictionary
35 
36  Output:
37  ofDict: OpenFoam dictionary
38 
39  Example:
40  We support two types of pyDict, one is to set the key value
41  as a list and have the value type as the first element.
42  This is needed in pyDAFoam.py to reminder users the date type
43  of a key value
44 
45  pyDict = {
46  "solverName": [str, "DASimpleFoam"],
47  "patches": [list, [1.0, 2.0]]
48  }
49 
50  The second type of pyDict is just like a usual dict:
51 
52  pyDict = {
53  "solverName": "DASimpleFoam",
54  "patches": [1.0, 2.0]
55  }
56 
57  Note: one can have multiple levels of sub dictionaries
58  */
59 
60  //PyObject_Print(pyDict,stdout,0);Info<<endl;
61 
62  // size of the pyOpions keys
63  Py_ssize_t dictSize = PyDict_Size(pyDict);
64  // all the keys
65  PyObject* keys = PyDict_Keys(pyDict);
66  // loop over all the keys in pyDict and assign their values
67  // to ofDict
68  for (label i = 0; i < dictSize; i++)
69  {
70  // the ith key
71  PyObject* keyI = PyList_GetItem(keys, i);
72  // convert it to UTF8 such that we can use it in C++
73  const char* keyUTF8 = PyUnicode_AsUTF8(keyI);
74 
75  //std::cout << "Key is "<<keyUTF8<<std::endl;
76  // the actual value of this key, NOTE: it is a list
77  // the 1st one is its type and the 2nd one is its value
78  // this is because we need to make sure users prescribe
79  // the right data type in pyDAFoam.py, so we set the
80  // value like this in pyDAFoam.py:
81  // defOptions = {
82  // "solverName": [str, "DASimpleFOam"] }
83  // however, the above has one exception when using subDict
84  // in this case, we only have one element
85  PyObject* value = PyDict_GetItem(pyDict, keyI);
86  const char* valueTypeTmp = Py_TYPE(value)->tp_name;
87  PyObject* value1;
88  if (word(valueTypeTmp) == "list")
89  {
90  // there are two possibilities for this case
91  // for nonSubDict case, we need to have pyDict format
92  // like this: { "solverName": [str, "DASimpleDAFoam"] }
93  // however, for subDicts, we dont specify the type property
94  PyObject* value0 = PyList_GetItem(value, 0);
95  const char* valueTypeTmp0 = Py_TYPE(value0)->tp_name;
96  if (word(valueTypeTmp0) == "type")
97  {
98  // this is nonSubdict case
99  // the second element of value is the actual value
100  value1 = PyList_GetItem(value, 1);
101  }
102  else
103  {
104  // it is for subdicts, so value1=value
105  value1 = value;
106  }
107  }
108  else
109  {
110  // if we only have one element, set value1 = value
111  value1 = value;
112  }
113 
114  //PyObject_Print(value1,stdout,0);Info<<endl;
115  // get the type of value1
116  const char* valueType = Py_TYPE(value1)->tp_name;
117  bool overwrite = true;
118  if (word(valueType) == "str")
119  {
120  const char* valSet = PyUnicode_AsUTF8(value1);
121  ofDict.add(keyUTF8, word(valSet), overwrite);
122  }
123  else if (word(valueType) == "int")
124  {
125  long valSet = PyLong_AsLong(value1);
126  ofDict.add(keyUTF8, label(valSet), overwrite);
127  }
128  else if (word(valueType) == "bool")
129  {
130  label valSet = PyObject_IsTrue(value1);
131  ofDict.add(keyUTF8, valSet, overwrite);
132  }
133  else if (word(valueType) == "float")
134  {
135  scalar valSet = PyFloat_AS_DOUBLE(value1);
136  ofDict.add(keyUTF8, valSet, overwrite);
137  }
138  else if (word(valueType) == "list")
139  {
140  // size of the list
141  Py_ssize_t listSize = PyList_Size(value1);
142 
143  //Info<<listSize<<endl;
144  // create OpenFOAM lists to hold this list
145  // we create all the possible types
146  scalarList valSetScalar;
147  valSetScalar.setSize(label(listSize));
148  labelList valSetLabel;
149  valSetLabel.setSize(label(listSize));
150  List<word> valSetWord;
151  valSetWord.setSize(label(listSize));
152 
153  // need to check what type of list this is
154  // by checking its first element
155  PyObject* tmp = PyList_GetItem(value1, 0);
156  const char* tmpType = Py_TYPE(tmp)->tp_name;
157  word tmpTypeWord = word(tmpType);
158  // assign value to the OpenFOAM list
159  for (label j = 0; j < listSize; j++)
160  {
161  PyObject* valueListJ = PyList_GetItem(value1, j);
162  if (tmpTypeWord == "str")
163  {
164  const char* valSet = PyUnicode_AsUTF8(valueListJ);
165  valSetWord[j] = word(valSet);
166  }
167  else if (tmpTypeWord == "int")
168  {
169  long valSet = PyLong_AsLong(valueListJ);
170  valSetLabel[j] = label(valSet);
171  }
172  else if (tmpTypeWord == "float")
173  {
174  scalar valSet = PyFloat_AS_DOUBLE(valueListJ);
175  valSetScalar[j] = valSet;
176  }
177  else if (tmpTypeWord == "bool")
178  {
179  label valSet = PyObject_IsTrue(valueListJ);
180  valSetLabel[j] = valSet;
181  }
182  else
183  {
184  FatalErrorIn("pyDict2OFDict") << "Type: <" << tmpTypeWord << "> for " << keyUTF8
185  << " list is not supported! Options are:"
186  << " str, int, bool, and float!"
187  << abort(FatalError);
188  }
189  }
190 
191  // add the list to the ofDict dict
192  if (tmpTypeWord == "str")
193  {
194  ofDict.add(keyUTF8, valSetWord, overwrite);
195  }
196  else if (tmpTypeWord == "int")
197  {
198  ofDict.add(keyUTF8, valSetLabel, overwrite);
199  }
200  else if (tmpTypeWord == "float")
201  {
202  ofDict.add(keyUTF8, valSetScalar, overwrite);
203  }
204  else if (tmpTypeWord == "bool")
205  {
206  ofDict.add(keyUTF8, valSetLabel, overwrite);
207  }
208  }
209  else if (word(valueType) == "dict")
210  {
211  // if its a subdict, recursely call this function
212  dictionary subDict;
213  DAUtility::pyDict2OFDict(value1, subDict);
214  ofDict.add(keyUTF8, subDict, overwrite);
215  }
216  else
217  {
218  FatalErrorIn("pyDict2OFDict") << "Type: " << valueType << " for " << keyUTF8
219  << " is not supported! Options are: str, int, float, bool, list, and dict!"
220  << abort(FatalError);
221  }
222  //std::cout << "My type is " << valueType << std::endl;
223  }
224 }
225 
227  Vec vecIn,
228  const word prefix)
229 {
230  /*
231  Description:
232  Read a vector in binary form
233 
234  Input:
235  vecIn: a Petsc vector to read values into (also output)
236  prefix: Name of the Petsc vector from disk
237 
238  Example:
239  If the vector storing in the disk reads: dFdWVector.bin
240  Then read the vector using:
241  Vec dFdW;
242  codes to initialize vector dFdW....
243  readVectorBinary(dFdW,"dFdwVector");
244  NOTE: the prefix does not include ".bin"
245  */
246 
247  std::ostringstream fileNameStream("");
248  fileNameStream << prefix << ".bin";
249  word fileName = fileNameStream.str();
250 
251  PetscViewer viewer;
252  PetscViewerBinaryOpen(PETSC_COMM_WORLD, fileName.c_str(), FILE_MODE_READ, &viewer);
253  VecLoad(vecIn, viewer);
254  PetscViewerDestroy(&viewer);
255 
256  return;
257 }
258 
260  const Vec vecIn,
261  const word prefix)
262 {
263  /*
264  Description:
265  Write a vector in binary form
266 
267  Input:
268  vecIn: a Petsc vector to write to disk
269  prefix: Name of the Petsc vector to write to disk
270 
271  Example:
272  Vec dFdW;
273  codes to initialize vector dFdW....
274  writeVectorBinary(dFdW,"dFdwVector");
275  This will write the dFdW vector to disk with name "dFdWVector.bin"
276  NOTE: the prefix does not include ".bin"
277  */
278 
279  std::ostringstream fileNameStream("");
280  fileNameStream << prefix << ".bin";
281  word fileName = fileNameStream.str();
282 
283  PetscViewer viewer;
284  PetscViewerBinaryOpen(PETSC_COMM_WORLD, fileName.c_str(), FILE_MODE_WRITE, &viewer);
285  VecView(vecIn, viewer);
286  PetscViewerDestroy(&viewer);
287 
288  return;
289 }
290 
292  const Vec vecIn,
293  const word prefix)
294 {
295  /*
296  Description:
297  Write a vector in ASCII form
298  Input:
299  vecIn: a Petsc vector to write to disk
300  prefix: Name of the Petsc vector to write to disk
301  Example:
302  Vec dFdW;
303  codes to initialize vector dFdW....
304  writeVectorASCII(dFdW,"dFdwVector");
305  This will write the dFdW vector to disk with name "dFdWVector.dat"
306  NOTE: the prefix does not include ".dat"
307  */
308 
309  std::ostringstream fileNameStream("");
310  fileNameStream << prefix << ".dat";
311  word fileName = fileNameStream.str();
312 
313  PetscViewer viewer;
314  PetscViewerASCIIOpen(PETSC_COMM_WORLD, fileName.c_str(), &viewer);
315  PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB); // write all the digits
316  VecView(vecIn, viewer);
317  PetscViewerDestroy(&viewer);
318 
319  return;
320 }
321 
323  Mat matIn,
324  const word prefix)
325 {
326  /*
327  Description:
328  Read a matrix in binary form
329 
330  Input:
331  matIn: a Petsc matrix to read values into (also output)
332  prefix: Name of the Petsc matrix from disk
333 
334  Example:
335  If the matrix storing in the disk reads: dRdWMat.bin
336  Then read the matrix using:
337  Mat dRdW;
338  codes to initialize matrix dRdW....
339  readMatrixBinary(dRdW,"dRdwVector");
340  NOTE: the prefix does not include ".bin"
341  */
342 
343  std::ostringstream fileNameStream("");
344  fileNameStream << prefix << ".bin";
345  word fileName = fileNameStream.str();
346 
347  PetscViewer viewer;
348  PetscViewerBinaryOpen(PETSC_COMM_WORLD, fileName.c_str(), FILE_MODE_READ, &viewer);
349  MatLoad(matIn, viewer);
350  PetscViewerDestroy(&viewer);
351 
352  return;
353 }
354 
356  const Mat matIn,
357  const word prefix)
358 {
359  /*
360  Description:
361  Write a matrix in binary form
362 
363  Input:
364  matIn: a Petsc matrix to write to disk
365  prefix: Name of the Petsc matrix to write to disk
366 
367  Example:
368  Mat dRdW;
369  codes to initialize matrix dRdW....
370  writeMatrixBinary(dRdW,"dRdWMat");
371  This will write the dRdW matrix to disk with name "dRdWMat.bin"
372  NOTE: the prefix does not include ".bin"
373  */
374 
375  std::ostringstream fileNameStream("");
376  fileNameStream << prefix << ".bin";
377  word fileName = fileNameStream.str();
378 
379  PetscViewer viewer;
380  PetscViewerBinaryOpen(PETSC_COMM_WORLD, fileName.c_str(), FILE_MODE_WRITE, &viewer);
381  MatView(matIn, viewer);
382  PetscViewerDestroy(&viewer);
383 
384  return;
385 }
386 
388  const Mat matIn,
389  const word prefix)
390 {
391  /*
392  Description:
393  Write a matrix in ASCII form
394 
395  Input:
396  matIn: a Petsc matrix to write to disk
397  prefix: Name of the Petsc matrix to write to disk
398 
399  Example:
400  Mat dRdW;
401  codes to initialize matrix dRdW....
402  writeMatrixASCII(dRdW,"dRdWMat");
403  This will write the dRdW matrix to disk with name "dRdWMat.dat"
404  NOTE: the prefix does not include ".dat"
405  */
406 
407  std::ostringstream fileNameStream("");
408  fileNameStream << prefix << ".dat";
409  word fileName = fileNameStream.str();
410 
411  PetscViewer viewer;
412  PetscViewerASCIIOpen(PETSC_COMM_WORLD, fileName.c_str(), &viewer);
413  MatView(matIn, viewer);
414  PetscViewerDestroy(&viewer);
415 
416  return;
417 }
418 
420  const dictionary& allOptions,
421  volScalarField& var,
422  const label printToScreen)
423 {
424  /*
425  Description:
426  Bound a field variable according to the bounds defined
427  in the allOptions dict. This is a overload function for volScalarField
428 
429  Input:
430  allOptions: a dictionary that has upper and lower bound values
431  We need to give specific name for the bounds, i.e.,
432 
433  variable name + Max -> setting upper bound
434 
435  variable name + Min -> setting lower bound
436 
437  Input & Output:
438  var: an OpenFOAM field variable to be bounded
439 
440  Example:
441  dictionary allOptions;
442  allOptions.set("pMax", 120000.0);
443  allOptions.set("pMin", 80000.0);
444  DAUtility daUtil;
445  volScalarField p = .... // initialize p
446  daUtil.boundVar(allOptions, p);
447  */
448 
449  const scalar vGreat = 1e200;
450  label useUpperBound = 0, useLowerBound = 0;
451 
452  dictionary varBoundsDict = allOptions.subDict("primalVarBounds");
453 
454  word lowerBoundName = var.name() + "Min";
455  word upperBoundName = var.name() + "Max";
456 
457  scalar varMin = varBoundsDict.lookupOrDefault<scalar>(lowerBoundName, -vGreat);
458  scalar varMax = varBoundsDict.lookupOrDefault<scalar>(upperBoundName, vGreat);
459 
460  forAll(var, cellI)
461  {
462  if (var[cellI] <= varMin)
463  {
464  var[cellI] = varMin;
465  useLowerBound = 1;
466  }
467  if (var[cellI] >= varMax)
468  {
469  var[cellI] = varMax;
470  useUpperBound = 1;
471  }
472  }
473 
474  forAll(var.boundaryField(), patchI)
475  {
476  forAll(var.boundaryField()[patchI], faceI)
477  {
478  if (var.boundaryFieldRef()[patchI][faceI] <= varMin)
479  {
480  var.boundaryFieldRef()[patchI][faceI] = varMin;
481  useLowerBound = 1;
482  }
483  if (var.boundaryFieldRef()[patchI][faceI] >= varMax)
484  {
485  var.boundaryFieldRef()[patchI][faceI] = varMax;
486  useUpperBound = 1;
487  }
488  }
489  }
490 
491  if (printToScreen)
492  {
493  if (useUpperBound)
494  {
495  Info << "Bounding " << var.name() << "<" << varMax << endl;
496  }
497  if (useLowerBound)
498  {
499  Info << "Bounding " << var.name() << ">" << varMin << endl;
500  }
501  }
502 
503  return;
504 }
505 
507  const dictionary& allOptions,
508  volVectorField& var,
509  const label printToScreen)
510 {
511  /*
512  Description:
513  Bound a field variable according to the bounds defined
514  in the allOptions dict. This is a overload function for volVectorField
515 
516  Input:
517  allOptions: a dictionary that has upper and lower bound values
518  We need to give specific name for the bounds, i.e.,
519 
520  variable name + Max -> setting upper bound
521 
522  variable name + Min -> setting lower bound
523 
524  Input & Output:
525  var: an OpenFOAM field variable to be bounded
526 
527  Example:
528  dictionary allOptions;
529  allOptions.set("pMax", 120000.0);
530  allOptions.set("pMin", 80000.0);
531  DAUtility daUtil;
532  volScalarField p = .... // initialize p
533  daUtil.boundVar(allOptions, p);
534  */
535 
536  const scalar vGreat = 1e200;
537  label useUpperBound = 0, useLowerBound = 0;
538 
539  dictionary varBoundsDict = allOptions.subDict("primalVarBounds");
540 
541  word lowerBoundName = var.name() + "Min";
542  word upperBoundName = var.name() + "Max";
543 
544  scalar varMin = varBoundsDict.lookupOrDefault<scalar>(lowerBoundName, -vGreat);
545  scalar varMax = varBoundsDict.lookupOrDefault<scalar>(upperBoundName, vGreat);
546 
547  forAll(var, cellI)
548  {
549  for (label i = 0; i < 3; i++)
550  {
551  if (var[cellI][i] <= varMin)
552  {
553  var[cellI][i] = varMin;
554  useLowerBound = 1;
555  }
556  if (var[cellI][i] >= varMax)
557  {
558  var[cellI][i] = varMax;
559  useUpperBound = 1;
560  }
561  }
562  }
563 
564  forAll(var.boundaryField(), patchI)
565  {
566  forAll(var.boundaryField()[patchI], faceI)
567  {
568  for (label i = 0; i < 3; i++)
569  {
570  if (var.boundaryFieldRef()[patchI][faceI][i] <= varMin)
571  {
572  var.boundaryFieldRef()[patchI][faceI][i] = varMin;
573  useLowerBound = 1;
574  }
575  if (var.boundaryFieldRef()[patchI][faceI][i] >= varMax)
576  {
577  var.boundaryFieldRef()[patchI][faceI][i] = varMax;
578  useUpperBound = 1;
579  }
580  }
581  }
582  }
583 
584  if (printToScreen)
585  {
586  if (useUpperBound)
587  {
588  Info << "Bounding " << var.name() << "<" << varMax << endl;
589  }
590  if (useLowerBound)
591  {
592  Info << "Bounding " << var.name() << ">" << varMin << endl;
593  }
594  }
595 
596  return;
597 }
598 
599 globalIndex DAUtility::genGlobalIndex(const label localIndexSize)
600 {
601  /*
602  Generate a glocal index system based on the local index size
603  such that we can use it to map a local index to a global one
604 
605  Input:
606  -----
607  localIndexSize: the SIZE of local index
608 
609  Output:
610  ------
611  globalIndex object: the global index object to map a local index
612  to a global index
613 
614  Example:
615  --------
616  If the local index reads:
617 
618  On processor 0:
619  labelList sampleList = {0, 1, 2};
620  globalIndex glbSample = genGlobalIndex(sampleList.size());
621 
622  On processor 1:
623  labelList sampleList = {0, 1};
624  globalIndex glbSample = genGlobalIndex(sampleList.size());
625 
626  After each processor calls genGlobalIndex and get the glbSample
627  object, we can use it to map a local index to a global one,
628  e.g., on processor 0, if we call:
629 
630  label glxIdx = glbSample.toGlobal(1);
631 
632  it will return glbIdx = 1;
633  However, on processor 1, if we call
634 
635  label glxIdx = glbSample.toGlobal(1);
636 
637  it will return glbIdx = 4;
638 
639  The date storage structure is illustrated as follows
640 
641  global index -> 0 1 2 3 4
642  local index -> 0 1 2 0 1
643  ----- ===
644  proc0 proc1
645 
646  */
647  globalIndex result(localIndexSize);
648  return result;
649 }
650 
652  const scalar val,
653  const scalar refVal,
654  const scalar tol)
655 {
656  /*
657  Description:
658  check whether a value is close to a reference value by a tolerance
659 
660  Input:
661  val: value to check
662  refVal: reference value
663  tol: tolerance
664 
665  Output:
666  return: 1 means fabs(val-refVal) < tol
667  */
668  if (fabs(val - refVal) < tol)
669  {
670  return 1;
671  }
672  else
673  {
674  return 0;
675  }
676 }
677 
678 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
679 
680 } // End namespace Foam
681 
682 // ************************************************************************* //
Foam::DAUtility::readMatrixBinary
static void readMatrixBinary(Mat matIn, const word prefix)
read petsc matrix in binary format
Definition: DAUtility.C:322
allOptions
const dictionary & allOptions
Definition: createRefsRhoSimpleC.H:15
Foam::DAUtility::readVectorBinary
static void readVectorBinary(Vec vecIn, const word prefix)
read petsc vector in binary format
Definition: DAUtility.C:226
Foam::DAUtility::DAUtility
DAUtility()
Constructors.
Definition: DAUtility.C:16
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DAUtility::~DAUtility
virtual ~DAUtility()
Destructor.
Definition: DAUtility.C:20
DAUtility.H
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:651
Foam::DAUtility::pyDict2OFDict
static void pyDict2OFDict(PyObject *pyDict, dictionary &ofDict)
convert a python dictionary object to OpenFoam dictionary
Definition: DAUtility.C:24
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:419
Foam::DAUtility::writeVectorBinary
static void writeVectorBinary(const Vec vecIn, const word prefix)
write petsc vector in binary format
Definition: DAUtility.C:259
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAUtility::writeMatrixBinary
static void writeMatrixBinary(const Mat matIn, const word prefix)
write petsc matrix in binary format
Definition: DAUtility.C:355
Foam::DAUtility::writeVectorASCII
static void writeVectorASCII(const Vec vecIn, const word prefix)
write petsc vector in ascii format
Definition: DAUtility.C:291
Foam::DAUtility::writeMatrixASCII
static void writeMatrixASCII(const Mat matIn, const word prefix)
write petsc matrix in ascii format
Definition: DAUtility.C:387
Foam::DAUtility::genGlobalIndex
static globalIndex genGlobalIndex(const label localIndexSize)
generate global index numbering for local-global index transferring
Definition: DAUtility.C:599