cyclicTransform.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cyclicTransform.H"
27 #include "unitConversion.H"
28 #include "IOmanip.H"
29 #include "stringOps.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 template<class Type>
37 Type sum(const Type& x, const bool global)
38 {
39  return global ? returnReduce(x, sumOp<Type>()) : x;
40 }
41 
42 template<class Type>
43 Type sum(const Field<Type>& x, const bool global)
44 {
45  return global ? gSum(x) : sum(x);
46 }
47 
48 template<class Type>
49 Type sum(const tmp<Field<Type>>& x, const bool global)
50 {
51  const Type s = sum(x(), global);
52  x.clear();
53  return s;
54 }
55 
56 }
57 
58 
59 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
64 }
65 
68 {
69  "unspecified",
70  "none",
71  "rotational",
72  "translational"
73 };
74 
76 {
77  "transformType",
78  "transform",
79  "rotationAxis",
80  "rotationCentre",
81  "rotationAngle",
82  "separation",
83  "separationVector"
84 };
85 
86 
87 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
88 
89 void Foam::cyclicTransform::update()
90 {
91  if (!transformComplete_)
92  {
93  return;
94  }
95 
96  switch (transformType_)
97  {
98  case UNSPECIFIED:
99  break;
100 
101  case NONE:
102  transform_ = transformer();
103  break;
104 
105  case ROTATIONAL:
106  if (rotationAngle_ == 0)
107  {
108  transform_ = transformer();
109  }
110  else
111  {
112  const tensor R =
113  quaternion(rotationAxis_, rotationAngle_).R();
114 
115  if (mag(rotationCentre_) == 0)
116  {
117  transform_ = transformer::rotation(R);
118  }
119  else
120  {
121  transform_ =
122  transformer::translation(rotationCentre_)
124  & transformer::translation(- rotationCentre_);
125  }
126  }
127  break;
128 
129  case TRANSLATIONAL:
130  if (mag(separation_) == 0)
131  {
132  transform_ = transformer();
133  }
134  else
135  {
136  transform_ = transformer::translation(separation_);
137  }
138  break;
139  }
140 }
141 
142 
143 bool Foam::cyclicTransform::set
144 (
145  const cyclicTransform& t,
146  const scalar lengthScale,
147  const scalar matchTolerance
148 )
149 {
150  // If the supplied transform is unspecified then there is nothing to do
151  if (t.transformType_ == UNSPECIFIED)
152  {
153  return true;
154  }
155 
156  // If this transform is specified then we need to check that it is
157  // compatible with the supplied transform
158  if (transformType_ != UNSPECIFIED)
159  {
160  // If the transforms are of different types then they are incompatible
161  if (transformType_ != t.transformType_)
162  {
163  return false;
164  }
165 
166  // Length-tolerance
167  const scalar lengthTolerance = lengthScale*matchTolerance;
168 
169  // If the transforms are both rotational then the axes must be in the
170  // same direction, the centre points must lie on the same line, and the
171  // angles (if available) must be opposing.
172  if (transformType_ == ROTATIONAL)
173  {
174  const scalar dot = rotationAxis_ & t.rotationAxis_;
175 
176  if (mag(dot) < 1 - matchTolerance)
177  {
178  return false;
179  }
180 
181  if
182  (
183  (rotationAxis_ & (rotationCentre_ - t.rotationCentre_))
184  > lengthTolerance
185  )
186  {
187  return false;
188  }
189 
190  if (transformComplete_ && t.transformComplete_)
191  {
192  if
193  (
194  mag(rotationAngle_ - sign(dot)*t.rotationAngle_)
195  > matchTolerance
196  )
197  {
198  return false;
199  }
200  }
201  }
202 
203  // If the transforms are both translational then the separations must
204  // be opposing
205  if (transformType_ == TRANSLATIONAL)
206  {
207  if (transformComplete_ && t.transformComplete_)
208  {
209  if (mag(separation_ - t.separation_) > lengthTolerance)
210  {
211  return false;
212  }
213  }
214  }
215  }
216 
217  // If the supplied transform is more complete then overwrite this with it
218  if (!transformComplete_ && t.transformComplete_)
219  {
220  *this = t;
221  }
222 
223  return true;
224 }
225 
226 
227 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
228 
230 :
231  cyclicTransform(true)
232 {}
233 
234 
236 (
237  const bool defaultIsNone
238 )
239 :
240  transformType_(defaultIsNone ? NONE : UNSPECIFIED),
241  rotationAxis_(vector::uniform(NaN)),
242  rotationCentre_(vector::uniform(NaN)),
243  rotationAngle_(NaN),
244  separation_(vector::uniform(NaN)),
245  transformComplete_(transformType_ == NONE),
246  transform_()
247 {}
248 
249 
251 (
252  const dictionary& dict,
253  const bool defaultIsNone
254 )
255 :
256  transformType_
257  (
258  transformTypeNames
259  [
260  dict.lookupOrDefaultBackwardsCompatible<word>
261  (
262  {"transformType", "transform"},
263  transformTypeNames[defaultIsNone ? NONE : UNSPECIFIED]
264  )
265  ]
266  ),
267  rotationAxis_
268  (
269  transformType_ == ROTATIONAL
270  ? normalised(dict.lookup<vector>("rotationAxis"))
271  : vector::uniform(NaN)
272  ),
273  rotationCentre_
274  (
275  transformType_ == ROTATIONAL
276  ? dict.lookup<point>("rotationCentre")
277  : point::uniform(NaN)
278  ),
279  rotationAngle_
280  (
281  dict.lookupOrDefault<scalar>("rotationAngle", unitDegrees, NaN)
282  ),
283  separation_
284  (
285  transformType_ == TRANSLATIONAL
286  ? (
288  (
289  {"separation", "separationVector"}
290  )
291  )
292  : vector::uniform(NaN)
293  ),
294  transformComplete_
295  (
296  transformType_ == NONE
297  || (
298  transformType_ == ROTATIONAL
299  && dict.found("rotationAngle")
300  )
301  || (
302  transformType_ == TRANSLATIONAL
303  && (dict.found("separation") || dict.found("separationVector"))
304  )
305  ),
306  transform_()
307 {
308  if (transformComplete_)
309  {
310  update();
311  }
312 }
313 
314 
316 (
317  const word& name,
318  const vectorField& areas,
319  const cyclicTransform& transform,
320  const word& nbrName,
321  const cyclicTransform& nbrTransform,
322  const scalar matchTolerance,
323  const bool global
324 )
325 :
327 {
328  // Calculate the total (vector) areas for the supplied patch data
329  const vector area = sum(areas, global);
330 
331  // Calculate patch length scales
332  const scalar lengthScale = sqrt(mag(area));
333 
334  // Copy as much data from the neighbour as possible
335  if (!transformComplete_ && nbrTransform.transformType_ != UNSPECIFIED)
336  {
337  if (!set(inv(nbrTransform), lengthScale, matchTolerance))
338  {
340  << "Patch " << name
341  << " and neighbour patch " << nbrName
342  << " have incompatible transforms:" << nl << nl << incrIndent;
343 
345  << indent << name << nl << indent << token::BEGIN_BLOCK << nl
346  << incrIndent;
347 
349 
351  << decrIndent << indent << token::END_BLOCK << nl << nl;
352 
354  << indent << nbrName << nl << indent << token::BEGIN_BLOCK << nl
355  << incrIndent;
356 
357  nbrTransform.write(FatalError);
358 
360  << decrIndent << indent << token::END_BLOCK << nl << nl;
361 
363  << decrIndent << exit(FatalError);
364  }
365  }
366 }
367 
368 
370 (
371  const word& name,
372  const pointField& ctrs,
373  const vectorField& areas,
374  const cyclicTransform& transform,
375  const word& nbrName,
376  const pointField& nbrCtrs,
377  const vectorField& nbrAreas,
378  const cyclicTransform& nbrTransform,
379  const scalar matchTolerance,
380  const bool global
381 )
382 :
384  (
385  name,
386  areas,
387  transform,
388  nbrName,
389  nbrTransform,
390  matchTolerance
391  )
392 {
393  // If there is no geometry from which to calculate the transform then
394  // nothing can be calculated
395  if (sum(areas.size(), global) == 0 || sum(nbrAreas.size(), global) == 0)
396  {
397  return;
398  }
399 
400  // Calculate the total (vector) areas for the supplied patch data
401  const vector area = sum(areas, global);
402  const vector nbrArea = sum(nbrAreas, global);
403 
404  // Calculate the centroids for the supplied patch data
405  const scalarField magAreas(mag(areas));
406  const scalarField magNbrAreas(mag(nbrAreas));
407  const scalar sumMagAreas = sum(magAreas, global);
408  const scalar sumMagNbrAreas = sum(magNbrAreas, global);
409  const point ctr = sum(ctrs*magAreas, global)/sumMagAreas;
410  const point nbrCtr = sum(nbrCtrs*magNbrAreas, global)/sumMagNbrAreas;
411 
412  // Calculate patch length scales
413  const scalar lengthScale = sqrt(sumMagAreas);
414 
415  // Calculate the transformation from the patch geometry
416  if (!transformComplete_)
417  {
418  // Store the old transformation type
419  const transformTypes oldTransformType = transformType_;
420 
421  // Calculate the average patch normals
422  const vector normal = normalised(area);
423  const vector negNbrNormal = - normalised(nbrArea);
424 
425  // Calculate the angle and distance separations
426  const scalar dot = normal & negNbrNormal;
427  const vector delta = ctr - nbrCtr;
428 
429  // Determine the type of transformation if it has not been specified
430  if (transformType_ == UNSPECIFIED)
431  {
432  transformType_ =
433  dot < 1 - rootSmall
434  ? ROTATIONAL
435  : mag(delta) > lengthScale*rootSmall
436  ? TRANSLATIONAL
437  : NONE;
438  }
439 
440  // If the transformation is known to be rotational, then we need to
441  // calculate the angle. If the transformation was previously
442  // unspecified then we also need to calculate the axis and the centre
443  // of rotation.
444  if (transformType_ == ROTATIONAL)
445  {
446  // Calculate the axis, if necessary
447  if (transformType_ != oldTransformType)
448  {
449  const vector midNormal = normalised(normal + negNbrNormal);
450  const vector axis =
451  (ctr - nbrCtr)
452  ^ (
453  normal*(negNbrNormal & midNormal)
454  - negNbrNormal*(normal & midNormal)
455  );
456  const vector axis180 =
457  (ctr - nbrCtr)
458  ^ (normal - negNbrNormal);
459 
460  rotationAxis_ =
461  normalised
462  (
463  mag(axis) > lengthScale*rootSmall
464  ? axis
465  : axis180
466  );
467  }
468 
469  const tensor PerpA = tensor::I - sqr(rotationAxis_);
470  const vector normalPerpA = normalised(PerpA & normal);
471  const vector negNbrNormalPerpA = normalised(PerpA & negNbrNormal);
472  const scalar theta =
473  acos(min(max(normalPerpA & negNbrNormalPerpA, -1), 1));
474  rotationAngle_ =
475  - sign((normalPerpA ^ negNbrNormalPerpA) & rotationAxis_)*theta;
476 
477  // Calculate the angle
478  // Calculate the centre of rotation, if necessary
479  if (transformType_ != oldTransformType)
480  {
481  const tensor R = quaternion(rotationAxis_, theta).R();
482  tensor A = tensor::I - R;
483  vector b = ctr - (R & nbrCtr);
484  const label i = findMax(cmptMag(rotationAxis_));
485  forAll(b, j)
486  {
487  A(i, j) = j == i;
488  }
489  b[i] = 0;
490  rotationCentre_ = inv(A) & b;
491  }
492  }
493 
494  // If the transformation is known to be translational then we just need
495  // to set the separation.
496  if (transformType_ == TRANSLATIONAL)
497  {
498  separation_ = delta;
499  }
500 
501  // Update the transform object
502  transformComplete_ = true;
503  update();
504 
505  // Print results of calculation
506  if (debug)
507  {
508  Info<< "Transformation calculated between patches " << name
509  << " and " << nbrName << ":" << nl << token::BEGIN_BLOCK << nl
510  << incrIndent;
511 
513 
515  }
516  }
517 
518  // Check the transformation is correct to within the matching tolerance
519  const point nbrCtrT =
520  transform_.transformPosition(nbrCtr);
521 
522  const scalar ctrNbrCtrTDistance = mag(ctr - nbrCtrT);
523 
524  if (ctrNbrCtrTDistance > lengthScale*matchTolerance)
525  {
527  str << "Patches " << name << " and " << nbrName << " are potentially "
528  << "not geometrically similar enough to be coupled." << nl << nl
529  << "The distance between the transformed centres of these patches "
530  << "is " << ctrNbrCtrTDistance << ", which is greater than the "
531  << "patch length scale (" << lengthScale << ") multiplied by the "
532  << "match tolerance (" << matchTolerance << ")." << nl << nl
533  << "Check that the patches are geometrically similar and that any "
534  << "transformations defined between them are correct." << nl << nl
535  << "If the patches and their transformations are defined correctly "
536  << "but small irregularities in the mesh mean this geometric test "
537  << "is failing, then it might be appropriate to relax the failure "
538  << "criteria by increasing the \"matchTolerance\" setting for "
539  << "these patches in the \"polyMesh/boundary\" file.";
541  << nl << stringOps::breakIntoIndentedLines(str.str()).c_str()
542  << exit(FatalError);
543  }
544 }
545 
546 
547 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
548 
550 {}
551 
552 
553 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
554 
556 {
557  const label oldPrecision = os.precision();
558 
559  os.precision(16);
560 
561  if (transformType_ != UNSPECIFIED)
562  {
563  writeEntry(os, "transformType", transformTypeNames[transformType_]);
564  }
565 
566  if (transformType_ == ROTATIONAL)
567  {
568  writeEntry(os, "rotationAxis", rotationAxis_);
569  writeEntry(os, "rotationCentre", rotationCentre_);
570 
571  if (transformComplete_)
572  {
573  writeEntry(os, "rotationAngle", unitDegrees, rotationAngle_);
574  }
575  }
576 
577  if (transformType_ == TRANSLATIONAL)
578  {
579  if (transformComplete_)
580  {
581  writeEntry(os, "separation", separation_);
582  }
583  }
584 
585  os.precision(oldPrecision);
586 }
587 
588 
590 {
591  OStringStream oss;
592 
593  switch (transformType_)
594  {
595  case UNSPECIFIED:
596  break;
597 
598  case NONE:
599  break;
600 
601  case ROTATIONAL:
602  {
603  oss << "axis=" << rotationAxis_
604  << ", centre=" << rotationCentre_
605  << ", angle=";
606 
607  if (transformComplete_)
608  {
609  oss << unitDegrees.toUser(rotationAngle_);
610  }
611  else
612  {
613  oss << '?';
614  }
615 
616  break;
617  }
618 
619  case TRANSLATIONAL:
620  {
621  oss << "separation=";
622 
623  if (transformComplete_)
624  {
625  oss << separation_;
626  }
627  else
628  {
629  oss << '?';
630  }
631 
632  break;
633  }
634  }
635 
636  return oss.str();
637 }
638 
639 
640 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
641 
642 Foam::cyclicTransform Foam::operator&
643 (
644  const transformer& t,
645  const cyclicTransform& c0
646 )
647 {
648  cyclicTransform c1(c0);
649 
650  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
651  {
652  c1.rotationAxis_ = normalised(t.transform(c1.rotationAxis_));
653  c1.rotationCentre_ = t.transformPosition(c1.rotationCentre_);
654  }
655 
656  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
657  {
658  if (c1.transformComplete_)
659  {
660  c1.separation_ = t.transform(c1.separation_);
661  }
662  }
663 
664  if (c1.transformComplete_)
665  {
666  c1.update();
667  }
668 
669  return c1;
670 }
671 
672 
674 {
675  cyclicTransform c1(c0);
676 
677  if (c1.transformType_ == cyclicTransform::ROTATIONAL)
678  {
679  if (c1.transformComplete_)
680  {
681  c1.rotationAngle_ = - c1.rotationAngle_;
682  }
683  }
684 
685  if (c1.transformType_ == cyclicTransform::TRANSLATIONAL)
686  {
687  if (c1.transformComplete_)
688  {
689  c1.separation_ = - c1.separation_;
690  }
691  }
692 
693  if (c1.transformComplete_)
694  {
695  c1.update();
696  }
697 
698  return c1;
699 }
700 
701 
702 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Pre-declare SubField and related Field type.
Definition: Field.H:83
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Output to memory buffer stream.
Definition: OStringStream.H:52
string str() const
Return the string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual int precision() const =0
Get precision of output field.
static const Tensor I
Definition: Tensor.H:83
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
Cyclic plane transformation.
string str() const
Generate a string representation of the transform.
static const NamedEnum< transformTypes, 4 > transformTypeNames
friend cyclicTransform inv(const cyclicTransform &c)
void write(Ostream &os) const
Write the data to a dictionary.
static const wordList keywords
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:751
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Quaternion class used to perform rotations in 3D space.
Definition: quaternion.H:61
tensor R() const
The rotation tensor corresponding the quaternion.
Definition: quaternionI.H:323
A class for handling character strings derived from std::string.
Definition: string.H:79
A class for managing temporary objects.
Definition: tmp.H:55
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
static transformer rotation(const tensor &T)
Construct a pure rotation transformer.
Definition: transformerI.H:43
static transformer translation(const vector &t)
Construct a pure translation transformer.
Definition: transformerI.H:31
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
T toUser(const T &) const
Convert a value to user units.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
static const coefficient A("A", dimPressure, 611.21)
string breakIntoIndentedLines(const string &str, const string::size_type nLength=80, const string::size_type nIndent=0)
Break a string up into indented lines.
Definition: stringOps.C:963
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void dot(LagrangianPatchField< typename innerProduct< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
Type gSum(const FieldField< Field, Type > &f)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
void cmptMag(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const unitConversion unitDegrees
dimensionedScalar acos(const dimensionedScalar &ds)
dictionary dict