GeometricScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "GeometricScalarField.H"
30 
31 #define TEMPLATE template<template<class> class PatchField, class GeoMesh>
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
41 template<template<class> class PatchField, class GeoMesh>
42 void stabilise
43 (
46  const dimensioned<scalar>& ds
47 )
48 {
50  stabilise(result.boundaryFieldRef(), gsf.boundaryField(), ds.value());
51 }
52 
53 
54 template<template<class> class PatchField, class GeoMesh>
55 tmp<GeometricField<scalar, PatchField, GeoMesh>> stabilise
56 (
57  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
58  const dimensioned<scalar>& ds
59 )
60 {
62  (
63  "stabilise(" + gsf.name() + ',' + ds.name() + ')',
64  gsf.mesh(),
65  ds.dimensions() + gsf.dimensions()
66  );
67 
68  stabilise(tRes.ref(), gsf, ds);
69 
70  return tRes;
71 }
72 
73 
74 template<template<class> class PatchField, class GeoMesh>
75 tmp<GeometricField<scalar, PatchField, GeoMesh>> stabilise
76 (
77  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
78  const dimensioned<scalar>& ds
79 )
80 {
81  const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
82 
83  tmp<GeometricField<scalar, PatchField, GeoMesh>> tRes
84  (
85  New
86  (
87  tgsf,
88  "stabilise(" + gsf.name() + ',' + ds.name() + ')',
89  ds.dimensions() + gsf.dimensions()
90  )
91  );
92 
93  stabilise(tRes.ref(), gsf, ds);
94 
95  tgsf.clear();
96 
97  return tRes;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
102 
103 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
104 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
105 
106 BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
107 BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
108 
109 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 template<template<class> class PatchField, class GeoMesh>
115 void pow
116 (
117  GeometricField<scalar, PatchField, GeoMesh>& Pow,
118  const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
119  const GeometricField<scalar, PatchField, GeoMesh>& gsf2
120 )
121 {
122  pow(Pow.primitiveFieldRef(), gsf1.primitiveField(), gsf2.primitiveField());
123  pow(Pow.boundaryFieldRef(), gsf1.boundaryField(), gsf2.boundaryField());
124 }
125 
126 
127 template<template<class> class PatchField, class GeoMesh>
128 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
129 (
130  const GeometricField<scalar, PatchField, GeoMesh>& f1,
131  const GeometricField<scalar, PatchField, GeoMesh>& f2
132 )
133 {
134  if
135  (
137  && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
138  )
139  {
141  << "pow() expects dimensionless parameters, but found" << nl;
142 
143  if (!f1.dimensions().dimensionless())
144  {
145  FatalError
146  << " Base field dimensions: " << f1.dimensions() << nl;
147  }
148  if (!f2.dimensions().dimensionless())
149  {
150  FatalError
151  << " Exponent field dimensions: " << f2.dimensions() << nl;
152  }
154  }
155 
157  (
158  "pow(" + f1.name() + ',' + f2.name() + ')',
159  f1.mesh(),
160  dimless
161  );
162 
163  pow(tresult.ref(), f1, f2);
165  return tresult;
166 }
167 
168 
169 template<template<class> class PatchField, class GeoMesh>
170 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
171 (
172  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
173  const GeometricField<scalar, PatchField, GeoMesh>& f2
174 )
175 {
176  const auto& f1 = tf1();
177 
178  if
179  (
181  && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
182  )
183  {
185  << "pow() expects dimensionless parameters, but found" << nl;
186 
187  if (!f1.dimensions().dimensionless())
188  {
189  FatalError
190  << " Base field dimensions: " << f1.dimensions() << nl;
191  }
192  if (!f2.dimensions().dimensionless())
193  {
194  FatalError
195  << " Exponent field dimensions: " << f2.dimensions() << nl;
196  }
198  }
199 
200  tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
201  (
202  New
203  (
204  tf1,
205  "pow(" + f1.name() + ',' + f2.name() + ')',
206  dimless
207  )
208  );
209 
210  pow(tresult.ref(), f1, f2);
211 
212  tf1.clear();
214  return tresult;
215 }
216 
217 
218 template<template<class> class PatchField, class GeoMesh>
219 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
220 (
221  const GeometricField<scalar, PatchField, GeoMesh>& f1,
222  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
223 )
224 {
225  const auto& f2 = tf2();
226 
227  if
228  (
230  && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
231  )
232  {
234  << "pow() expects dimensionless parameters, but found" << nl;
235 
236  if (!f1.dimensions().dimensionless())
237  {
238  FatalError
239  << " Base field dimensions: " << f1.dimensions() << nl;
240  }
241  if (!f2.dimensions().dimensionless())
242  {
243  FatalError
244  << " Exponent field dimensions: " << f2.dimensions() << nl;
245  }
247  }
248 
249  tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
250  (
251  New
252  (
253  tf2,
254  "pow(" + f1.name() + ',' + f2.name() + ')',
255  dimless
256  )
257  );
258 
259  pow(tresult.ref(), f1, f2);
260 
261  tf2.clear();
262 
263  return tresult;
264 }
265 
266 template<template<class> class PatchField, class GeoMesh>
267 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
268 (
269  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
270  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
271 )
272 {
273  const auto& f1 = tf1();
274  const auto& f2 = tf2();
275 
276  if
277  (
279  && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
280  )
281  {
283  << "pow() expects dimensionless parameters, but found" << nl;
284 
285  if (!f1.dimensions().dimensionless())
286  {
287  FatalError
288  << " Base field dimensions: " << f1.dimensions() << nl;
289  }
290  if (!f2.dimensions().dimensionless())
291  {
292  FatalError
293  << " Exponent field dimensions: " << f2.dimensions() << nl;
294  }
296  }
297 
298  tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
299  (
300  reuseTmpTmpGeometricField
301  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
302  (
303  tf1,
304  tf2,
305  "pow(" + f1.name() + ',' + f2.name() + ')',
306  dimless
307  )
308  );
309 
310  pow(tresult.ref(), f1, f2);
311 
312  tf1.clear();
313  tf2.clear();
315  return tresult;
316 }
317 
318 
319 template<template<class> class PatchField, class GeoMesh>
320 void pow
321 (
322  GeometricField<scalar, PatchField, GeoMesh>& tPow,
323  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
324  const dimensioned<scalar>& ds
325 )
326 {
327  pow(tPow.primitiveFieldRef(), gsf.primitiveField(), ds.value());
328  pow(tPow.boundaryFieldRef(), gsf.boundaryField(), ds.value());
329 }
330 
331 
332 template<template<class> class PatchField, class GeoMesh>
333 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
334 (
335  const GeometricField<scalar, PatchField, GeoMesh>& f1,
336  const dimensionedScalar& ds
337 )
338 {
339  if
340  (
342  && (!ds.dimensions().dimensionless())
343  )
344  {
346  << "pow() expects dimensionless parameters, but found" << nl
347  << " Exponent dimensions: " << ds.dimensions() << nl
348  << exit(FatalError);
349  }
350 
352  (
353  "pow(" + f1.name() + ',' + ds.name() + ')',
354  f1.mesh(),
355  pow(f1.dimensions(), ds)
356  );
357 
358  pow(tresult.ref(), f1, ds);
359 
360  return tresult;
361 }
362 
363 template<template<class> class PatchField, class GeoMesh>
364 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
365 (
366  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
367  const dimensionedScalar& ds
368 )
369 {
370  const auto& f1 = tf1();
371 
372  if
373  (
375  && (!ds.dimensions().dimensionless())
376  )
377  {
379  << "pow() expects dimensionless parameters, but found" << nl
380  << " Exponent dimensions: " << ds.dimensions() << nl
381  << exit(FatalError);
382  }
383 
384  tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
385  (
386  New
387  (
388  tf1,
389  "pow(" + f1.name() + ',' + ds.name() + ')',
390  pow(f1.dimensions(), ds)
391  )
392  );
393 
394  pow(tresult.ref(), f1, ds);
395 
396  tf1.clear();
397 
398  return tresult;
399 }
400 
401 template<template<class> class PatchField, class GeoMesh>
402 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
403 (
404  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
405  const scalar& s
406 )
407 {
408  return pow(gsf, dimensionedScalar(s));
409 }
410 
411 template<template<class> class PatchField, class GeoMesh>
412 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
413 (
414  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
415  const scalar& s
416 )
417 {
418  return pow(tgsf, dimensionedScalar(s));
419 }
420 
421 
422 template<template<class> class PatchField, class GeoMesh>
423 void pow
424 (
425  GeometricField<scalar, PatchField, GeoMesh>& tPow,
426  const dimensioned<scalar>& ds,
427  const GeometricField<scalar, PatchField, GeoMesh>& gsf
428 )
429 {
430  pow(tPow.primitiveFieldRef(), ds.value(), gsf.primitiveField());
431  pow(tPow.boundaryFieldRef(), ds.value(), gsf.boundaryField());
432 }
433 
434 
435 template<template<class> class PatchField, class GeoMesh>
436 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
437 (
438  const dimensionedScalar& ds,
439  const GeometricField<scalar, PatchField, GeoMesh>& f2
440 )
441 {
442  if
443  (
445  && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
446  )
447  {
449  << "pow() expects dimensionless parameters, but found" << nl;
450 
451  if (!ds.dimensions().dimensionless())
452  {
453  FatalError
454  << " Base scalar dimensions: " << ds.dimensions() << nl;
455  }
456  if (!f2.dimensions().dimensionless())
457  {
458  FatalError
459  << " Exponent field dimensions: " << f2.dimensions() << nl;
460  }
462  }
463 
465  (
466  "pow(" + ds.name() + ',' + f2.name() + ')',
467  f2.mesh(),
468  dimless
469  );
470 
471  pow(tresult.ref(), ds, f2);
473  return tresult;
474 }
475 
476 
477 template<template<class> class PatchField, class GeoMesh>
478 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
479 (
480  const dimensionedScalar& ds,
481  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
482 )
483 {
484  const auto& f2 = tf2();
485 
486  if
487  (
489  && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
490  )
491  {
493  << "pow() expects dimensionless parameters, but found" << nl;
494 
495  if (!ds.dimensions().dimensionless())
496  {
497  FatalError
498  << " Base scalar dimensions: " << ds.dimensions() << nl;
499  }
500  if (!f2.dimensions().dimensionless())
501  {
502  FatalError
503  << " Exponent field dimensions: " << f2.dimensions() << nl;
504  }
506  }
507 
508  tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
509  (
510  New
511  (
512  tf2,
513  "pow(" + ds.name() + ',' + f2.name() + ')',
514  dimless
515  )
516  );
517 
518  pow(tresult.ref(), ds, f2);
519 
520  tf2.clear();
521 
522  return tresult;
523 }
524 
525 template<template<class> class PatchField, class GeoMesh>
526 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
527 (
528  const scalar& s,
529  const GeometricField<scalar, PatchField, GeoMesh>& gsf
530 )
531 {
532  return pow(dimensionedScalar(s), gsf);
533 }
534 
535 template<template<class> class PatchField, class GeoMesh>
536 tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
537 (
538  const scalar& s,
539  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
540 )
541 {
542  return pow(dimensionedScalar(s), tgsf);
543 }
544 
545 
546 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
547 
548 template<template<class> class PatchField, class GeoMesh>
549 void atan2
550 (
551  GeometricField<scalar, PatchField, GeoMesh>& Atan2,
552  const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
553  const GeometricField<scalar, PatchField, GeoMesh>& gsf2
554 )
555 {
556  atan2
557  (
558  Atan2.primitiveFieldRef(),
559  gsf1.primitiveField(),
560  gsf2.primitiveField()
561  );
562  atan2
563  (
564  Atan2.boundaryFieldRef(),
565  gsf1.boundaryField(),
566  gsf2.boundaryField()
567  );
568 }
569 
570 
571 template<template<class> class PatchField, class GeoMesh>
572 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
573 (
574  const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
575  const GeometricField<scalar, PatchField, GeoMesh>& gsf2
576 )
577 {
579  (
580  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
581  gsf1.mesh(),
582  atan2(gsf1.dimensions(), gsf2.dimensions())
583  );
584 
585  atan2(tAtan2.ref(), gsf1, gsf2);
587  return tAtan2;
588 }
589 
590 
591 template<template<class> class PatchField, class GeoMesh>
592 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
593 (
594  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
595  const GeometricField<scalar, PatchField, GeoMesh>& gsf2
596 )
597 {
598  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
599 
600  tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
601  (
602  New
603  (
604  tgsf1,
605  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
606  atan2(gsf1.dimensions(), gsf2.dimensions())
607  )
608  );
609 
610  atan2(tAtan2.ref(), gsf1, gsf2);
611 
612  tgsf1.clear();
614  return tAtan2;
615 }
616 
617 
618 template<template<class> class PatchField, class GeoMesh>
619 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
620 (
621  const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
622  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
623 )
624 {
625  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
626 
627  tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
628  (
629  New
630  (
631  tgsf2,
632  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
633  atan2( gsf1.dimensions(), gsf2.dimensions())
634  )
635  );
636 
637  atan2(tAtan2.ref(), gsf1, gsf2);
638 
639  tgsf2.clear();
640 
641  return tAtan2;
642 }
643 
644 template<template<class> class PatchField, class GeoMesh>
645 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
646 (
647  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
648  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
649 )
650 {
651  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
652  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
653 
654  tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
655  (
656  reuseTmpTmpGeometricField
657  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
658  (
659  tgsf1,
660  tgsf2,
661  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
662  atan2(gsf1.dimensions(), gsf2.dimensions())
663  )
664  );
665 
666  atan2(tAtan2.ref(), gsf1, gsf2);
667 
668  tgsf1.clear();
669  tgsf2.clear();
671  return tAtan2;
672 }
673 
674 
675 template<template<class> class PatchField, class GeoMesh>
676 void atan2
677 (
678  GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
679  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
680  const dimensioned<scalar>& ds
681 )
682 {
683  atan2(tAtan2.primitiveFieldRef(), gsf.primitiveField(), ds.value());
684  atan2(tAtan2.boundaryFieldRef(), gsf.boundaryField(), ds.value());
685 }
686 
687 
688 template<template<class> class PatchField, class GeoMesh>
689 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
690 (
691  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
692  const dimensionedScalar& ds
693 )
694 {
696  (
697  "atan2(" + gsf.name() + ',' + ds.name() + ')',
698  gsf.mesh(),
699  atan2(gsf.dimensions(), ds)
700  );
701 
702  atan2(tAtan2.ref(), gsf, ds);
703 
704  return tAtan2;
705 }
706 
707 template<template<class> class PatchField, class GeoMesh>
708 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
709 (
710  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
711  const dimensionedScalar& ds
712 )
713 {
714  const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
715 
716  tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
717  (
718  New
719  (
720  tgsf,
721  "atan2(" + gsf.name() + ',' + ds.name() + ')',
722  atan2(gsf.dimensions(), ds)
723  )
724  );
725 
726  atan2(tAtan2.ref(), gsf, ds);
727 
728  tgsf.clear();
729 
730  return tAtan2;
731 }
732 
733 template<template<class> class PatchField, class GeoMesh>
734 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
735 (
736  const GeometricField<scalar, PatchField, GeoMesh>& gsf,
737  const scalar& s
738 )
739 {
740  return atan2(gsf, dimensionedScalar(s));
741 }
742 
743 template<template<class> class PatchField, class GeoMesh>
744 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
745 (
746  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
747  const scalar& s
748 )
749 {
750  return atan2(tgsf, dimensionedScalar(s));
751 }
752 
753 
754 template<template<class> class PatchField, class GeoMesh>
755 void atan2
756 (
757  GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
758  const dimensioned<scalar>& ds,
759  const GeometricField<scalar, PatchField, GeoMesh>& gsf
760 )
761 {
762  atan2(tAtan2.primitiveFieldRef(), ds.value(), gsf.primitiveField());
763  atan2(tAtan2.boundaryFieldRef(), ds.value(), gsf.boundaryField());
764 }
765 
766 
767 template<template<class> class PatchField, class GeoMesh>
768 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
769 (
770  const dimensionedScalar& ds,
771  const GeometricField<scalar, PatchField, GeoMesh>& gsf
772 )
773 {
775  (
776  "atan2(" + ds.name() + ',' + gsf.name() + ')',
777  gsf.mesh(),
778  atan2(ds, gsf.dimensions())
779  );
780 
781  atan2(tAtan2.ref(), ds, gsf);
783  return tAtan2;
784 }
785 
786 
787 template<template<class> class PatchField, class GeoMesh>
788 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
789 (
790  const dimensionedScalar& ds,
791  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
792 )
793 {
794  const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
795 
796  tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
797  (
798  New
799  (
800  tgsf,
801  "atan2(" + ds.name() + ',' + gsf.name() + ')',
802  atan2(ds, gsf.dimensions())
803  )
804  );
805 
806  atan2(tAtan2.ref(), ds, gsf);
807 
808  tgsf.clear();
809 
810  return tAtan2;
811 }
812 
813 template<template<class> class PatchField, class GeoMesh>
814 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
815 (
816  const scalar& s,
817  const GeometricField<scalar, PatchField, GeoMesh>& gsf
818 )
819 {
820  return atan2(dimensionedScalar(s), gsf);
821 }
822 
823 template<template<class> class PatchField, class GeoMesh>
824 tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
825 (
826  const scalar& s,
827  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
828 )
829 {
830  return atan2(dimensionedScalar(s), tgsf);
831 }
832 
833 
834 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
835 
836 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
837 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
838 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
839 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
840 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
841 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
842 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
843 UNARY_FUNCTION(scalar, scalar, sign, sign)
844 UNARY_FUNCTION(scalar, scalar, pos, pos)
845 UNARY_FUNCTION(scalar, scalar, pos0, pos0)
846 UNARY_FUNCTION(scalar, scalar, neg, neg)
847 UNARY_FUNCTION(scalar, scalar, neg0, neg0)
848 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
849 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
850 
851 UNARY_FUNCTION(scalar, scalar, exp, trans)
852 UNARY_FUNCTION(scalar, scalar, log, trans)
853 UNARY_FUNCTION(scalar, scalar, log10, trans)
854 UNARY_FUNCTION(scalar, scalar, sin, trans)
855 UNARY_FUNCTION(scalar, scalar, cos, trans)
856 UNARY_FUNCTION(scalar, scalar, tan, trans)
857 UNARY_FUNCTION(scalar, scalar, asin, trans)
858 UNARY_FUNCTION(scalar, scalar, acos, trans)
859 UNARY_FUNCTION(scalar, scalar, atan, trans)
860 UNARY_FUNCTION(scalar, scalar, sinh, trans)
861 UNARY_FUNCTION(scalar, scalar, cosh, trans)
862 UNARY_FUNCTION(scalar, scalar, tanh, trans)
863 UNARY_FUNCTION(scalar, scalar, asinh, trans)
864 UNARY_FUNCTION(scalar, scalar, acosh, trans)
865 UNARY_FUNCTION(scalar, scalar, atanh, trans)
866 UNARY_FUNCTION(scalar, scalar, erf, trans)
867 UNARY_FUNCTION(scalar, scalar, erfc, trans)
868 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
869 UNARY_FUNCTION(scalar, scalar, j0, trans)
870 UNARY_FUNCTION(scalar, scalar, j1, trans)
871 UNARY_FUNCTION(scalar, scalar, y0, trans)
872 UNARY_FUNCTION(scalar, scalar, y1, trans)
873 
874 
875 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
876 
877 #define BesselFunc(func) \
878  \
879 template<template<class> class PatchField, class GeoMesh> \
880 void func \
881 ( \
882  GeometricField<scalar, PatchField, GeoMesh>& gsf, \
883  const int n, \
884  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 \
885 ) \
886 { \
887  func(gsf.primitiveFieldRef(), n, gsf1.primitiveField()); \
888  func(gsf.boundaryFieldRef(), n, gsf1.boundaryField()); \
889 } \
890  \
891 template<template<class> class PatchField, class GeoMesh> \
892 tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
893 ( \
894  const int n, \
895  const GeometricField<scalar, PatchField, GeoMesh>& gsf \
896 ) \
897 { \
898  if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \
899  { \
900  FatalErrorInFunction \
901  << "Field is not dimensionless: " << gsf.dimensions() << nl \
902  << abort(FatalError); \
903  } \
904  \
905  auto tFunc = GeometricField<scalar, PatchField, GeoMesh>::New \
906  ( \
907  #func "(" + gsf.name() + ')', \
908  gsf.mesh(), \
909  dimless \
910  ); \
911  \
912  func(tFunc.ref(), n, gsf); \
913  \
914  return tFunc; \
915 } \
916  \
917 template<template<class> class PatchField, class GeoMesh> \
918 tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
919 ( \
920  const int n, \
921  const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf \
922 ) \
923 { \
924  const auto& gsf = tgsf(); \
925  \
926  if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \
927  { \
928  FatalErrorInFunction \
929  << "Field is not dimensionless: " << gsf.dimensions() << nl \
930  << abort(FatalError); \
931  } \
932  \
933  tmp<GeometricField<scalar, PatchField, GeoMesh>> tFunc \
934  ( \
935  New \
936  ( \
937  tgsf, \
938  #func "(" + gsf.name() + ')', \
939  dimless \
940  ) \
941  ); \
942  \
943  func(tFunc.ref(), n, gsf); \
944  \
945  tgsf.clear(); \
946  \
947  return tFunc; \
948 }
949 
950 BesselFunc(jn)
951 BesselFunc(yn)
952 
953 #undef BesselFunc
954 
955 
956 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
957 
958 } // End namespace Foam
959 
960 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
961 
962 #include "undefFieldFunctionsM.H"
963 
964 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
const Type & value() const noexcept
Return const reference to value.
dimensionedScalar acos(const dimensionedScalar &ds)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y0(const dimensionedScalar &ds)
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
const dimensionSet dimless
Dimensionless.
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
Scalar specific part of the implementation of GeometricField.
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
#define BesselFunc(func)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
Definition: dimensionSet.C:471
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
Definition: dimensionSet.H:241
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar atan(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
dimensionedScalar pow6(const dimensionedScalar &ds)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
dimensionedScalar cosh(const dimensionedScalar &ds)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)