60inline static void add_blockage_c
69 if (
b >
pars.blockageNoCT)
83inline static void add_blockage_f
93 a = 1.0 - (1.0 - a) * (1.0 -
b);
95 else if (
b <
pars.blockedFacePar || isHole)
122 <<
"No PDRblock set" <<
nl
126 const PDRblock& pdrBlock =
block();
127 const PDRblock::location& xgrid = pdrBlock.
grid().x();
128 const PDRblock::location& ygrid = pdrBlock.
grid().y();
129 const PDRblock::location& zgrid = pdrBlock.
grid().z();
135 int cxmin, cxmax, cymin, cymax, czmin, czmax;
136 int cfxmin, cfxmax, cfymin, cfymax, cfzmin, cfzmax;
143 rad_a = rad_b = 0.5*obs.
dia();
158 obs.
x() - rad_a, obs.
x() + rad_a,
160 xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
161 ); assert(cxmax >=0);
165 obs.
y() - rad_b, obs.
y() + rad_b,
167 yoverlap, &cymin, &cymax, &cfymin, &cfymax
168 ); assert(cymax >=0);
172 obs.
z(), obs.
z() + obs.
len(),
174 zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
175 ); assert(czmax >=0);
181 obs.
x(), obs.
y(), obs.
dia(),
183 pdrBlock.
grid().x(), cxmin, cfxmax,
184 pdrBlock.
grid().y(), cymin, cfymax,
190 for (label ix = cxmin; ix <= cfxmax; ix++)
192 for (label iy = cymin; iy <= cfymax; iy++)
194 for (label iz = czmin; iz <= cfzmax; iz++)
196 const scalar vol_block =
aboverlap(ix,iy) * zoverlap[iz];
197 v_block(ix,iy,iz) += vol_block;
204 if (!
pars.
two_d && (iz == czmin || iz == czmax))
207 const scalar both_ends_fac = (czmin == czmax ? 2.0 : 1.0);
210 * xgrid.
width(ix) * ygrid.
width(iy) * both_ends_fac;
213 const scalar temp =
c_count(ix,iy) * zoverlap[iz];
218 if (!
pars.
two_d && (iz == czmin || iz == czmax))
221 const scalar both_ends_fac = (czmin == czmax ? 2.0 : 1.0);
223 sub_count(ix,iy,iz).z() -= temp * both_ends_fac / 2.0;
233 area_block.x() =
ac_lblock(ix,iy) * zoverlap[iz];
234 area_block.y() =
bc_lblock(ix,iy) * zoverlap[iz];
293 for (label iz = czmin+1; iz < czmax; ++iz)
319 obs.
x(), obs.
x() + obs.
len(),
321 xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
322 ); assert(cxmax >=0);
326 obs.
y() - rad_a, obs.
y() + rad_a,
328 yoverlap, &cymin, &cymax, &cfymin, &cfymax
329 ); assert(cymax >=0);
333 obs.
z() - rad_b, obs.
z() + rad_b,
335 zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
336 ); assert(czmax >=0);
340 obs.
y(), obs.
z(), obs.
dia(),
342 pdrBlock.
grid().y(), cymin, cfymax,
343 pdrBlock.
grid().z(), czmin, cfzmax,
349 for (label iy = cymin; iy <= cfymax; iy++)
351 for (label iz = czmin; iz <= cfzmax; iz++)
353 for (label ix = cxmin; ix <= cxmax; ix++)
355 const scalar vol_block =
aboverlap(iy,iz) * xoverlap[ix];
356 v_block(ix,iy,iz) += vol_block;
359 if (ix == cxmin || ix == cxmax)
362 const scalar both_ends_fac = (cxmin == cxmax ? 2.0 : 1.0);
365 * ygrid.
width(iy) * zgrid.
width(iz) * both_ends_fac;
368 const scalar temp =
c_count(iy,iz) * xoverlap[ix];
372 if (ix == cfxmin || ix == cfxmax)
375 const scalar both_ends_fac = (cfxmin == cfxmax ? 2.0 : 1.0);
377 sub_count(ix,iy,iz).x() -= temp * both_ends_fac / 2.0;
380 area_block.y() =
ac_lblock(iy,iz) * xoverlap[ix];
381 area_block.z() =
bc_lblock(iy,iz) * xoverlap[ix];
430 for (label ix = cxmin+1; ix < cxmax; ++ix)
451 obs.
x() - rad_b, obs.
x() + rad_b,
453 xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
454 ); assert(cxmax >=0);
458 obs.
y(), obs.
y() + obs.
len(),
460 yoverlap, &cymin, &cymax, &cfymin, &cfymax
461 ); assert(cymax >=0);
465 obs.
z() - rad_a, obs.
z() + rad_a,
467 zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
468 ); assert(czmax >=0);
472 obs.
z(), obs.
x(), obs.
dia(),
474 pdrBlock.
grid().z(), czmin, cfzmax,
475 pdrBlock.
grid().x(), cxmin, cfxmax,
481 for (label iz = czmin; iz <= cfzmax; iz++)
483 for (label ix = cxmin; ix <= cfxmax; ix++)
485 for (label iy = cymin; iy <= cymax; iy++)
487 const scalar vol_block =
aboverlap(iz,ix) * yoverlap[iy];
488 v_block(ix,iy,iz) += vol_block;
491 if (iy == cymin || iy == cymax)
494 const scalar both_ends_fac = (cymin == cymax ? 2.0 : 1.0);
497 * zgrid.
width(iz) * xgrid.
width(ix) * both_ends_fac;
500 const scalar temp =
c_count(iz,ix) * yoverlap[iy];
505 if (iy == cfymin || iy == cfymax)
508 const scalar both_ends_fac = (cfymin == cfymax ? 2.0 : 1.0);
510 sub_count(ix,iy,iz).y() -= temp * both_ends_fac / 2.0;
513 area_block.z() =
ac_lblock(iz,ix) * yoverlap[iy];
514 area_block.x() =
bc_lblock(iz,ix) * yoverlap[iy];
563 for (label iy = cymin+1; iy < cymax; ++iy)
582 Info<<
"Unexpected orientation " << int(obs.
orient) <<
nl;
601 || (volumeSign < 0 && obs.vbkge >= 0)
602 || (volumeSign > 0 && obs.
vbkge < 0)
611 <<
"No PDRblock set" <<
nl
646 const auto spc = identifier.find_first_of(
" \t\n\v\f\r");
653 [&](
const PDRpatchDef&
p){
return patchName ==
p.patchName; },
678 if (obs.
span.
x() < 0.01)
682 else if (obs.
span.
y() < 0.01)
686 else if (obs.
span.
z() < 0.01)
693 <<
"Blowoff panel should have zero thickness" <<
nl
699 int cxmin, cxmax, cfxmin, cfxmax;
702 obs.
x(), obs.
x() + obs.
span.
x(),
704 xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
705 ); assert(cxmax >=0);
707 int cymin, cymax, cfymin, cfymax;
710 obs.
y(), obs.
y() + obs.
span.
y(),
712 yoverlap, &cymin, &cymax, &cfymin, &cfymax
713 ); assert(cymax >=0);
715 int czmin, czmax, cfzmin, cfzmax;
718 obs.
z(), obs.
z() + obs.
span.
z(),
720 zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
721 ); assert(czmax >=0);
723 two_d_overlap(xoverlap, cxmin, cxmax, yoverlap, cymin, cymax, aboverlap);
726 const scalar vbkge = obs.
vbkge;
727 const scalar xbkge = obs.
xbkge;
728 const scalar ybkge = obs.
ybkge;
729 const scalar zbkge = obs.
zbkge;
734 ((cxmin == cxmax) ? 0.5 : 1.0),
735 ((cymin == cymax) ? 0.5 : 1.0),
736 ((czmin == czmax) ? 0.5 : 1.0)
739 for (label ix = cxmin; ix <= cfxmax; ix++)
741 const scalar xov = xoverlap[ix];
743 scalar
area, cell_area, temp;
745 for (label iy = cymin; iy <= cfymax; iy++)
747 const scalar yov = yoverlap[iy];
749 for (label iz = czmin; iz <= cfzmax; iz++)
751 const scalar zov = zoverlap[iz];
758 (indir == -1 && ix == cfxmin)
759 || (indir == 1 && ix == cfxmax)
760 || (indir == -2 && iy == cfymin)
761 || (indir == 2 && iy == cfymax)
762 || (indir == -3 && iz == cfzmin)
763 || (indir == 3 && iz == cfzmax)
780 if (yov * zov >
pars.blockedFacePar)
782 face_patch(ix,iy,iz).x() = patchNum;
789 if (zov * xov >
pars.blockedFacePar)
791 face_patch(ix,iy,iz).y() = patchNum;
798 if (xov * yov >
pars.blockedFacePar)
800 face_patch(ix,iy,iz).z() = patchNum;
807 const scalar vol_block = aboverlap(ix,iy) * zov * vbkge;
808 v_block(ix,iy,iz) += vol_block;
811 if ((ix > cxmin && ix <= cxmax) || (cxmin == cxmax && ix == cfxmax))
813 temp = yov * zov * xbkge;
819 hole_in_face(ix,iy,iz).x() =
true;
821 add_blockage_f(face_block(ix,iy,iz).
x(), temp, hole_in_face(ix,iy,iz).
x());
822 if (temp >
pars.blockedFacePar && ! hole_in_face(ix,iy,iz).
x() && !isPatch)
828 if ((iy > cymin && iy <= cymax) || (cymin == cymax && iy == cfymax))
830 temp = zov * xov * ybkge;
833 hole_in_face(ix,iy,iz).y() =
true;
835 add_blockage_f(face_block(ix,iy,iz).
y(), temp, hole_in_face(ix,iy,iz).
y());
836 if (temp >
pars.blockedFacePar && ! hole_in_face(ix,iy,iz).
y() && !isPatch)
841 if ((iz > czmin && iz <= czmax) || (czmin == czmax && iz == cfzmax))
843 temp = xov * yov * zbkge;
846 hole_in_face(ix,iy,iz).z() =
true;
848 add_blockage_f(face_block(ix,iy,iz).z(), temp, hole_in_face(ix,iy,iz).z());
849 if (temp >
pars.blockedFacePar && ! hole_in_face(ix,iy,iz).z() && !isPatch)
858 area = yov * zov * xbkge;
859 if (ix < cxmin || ix > cxmax)
861 else if (ix > cxmin && ix < cxmax && xbkge >= 1.0)
863 along_block(ix,iy,iz).x() +=
area;
864 betai_inv1(ix,iy,iz).x() += vol_block / (1.0 -
area +
floatSMALL);
866 else if (ix == cxmin || ix == cxmax)
870 const scalar double_f = (cxmin == cxmax ? 1.0 : 0.5);
872 add_blockage_c(area_block_s(ix,iy,iz).
x(), dirn_block(ix,iy,iz).
x(), area, double_f);
873 cell_area = (ygrid.
width(iy) * zgrid.
width(iz));
874 surf(ix,iy,iz) += double_f *
area * cell_area;
875 betai_inv1(ix,iy,iz).x() += vol_block / (1.0 -
area +
floatSMALL);
882 if (temp > 0.0) { grating_count(ix,iy,iz).x() += temp; }
888 area = zov * xov * ybkge;
889 if (iy < cymin || iy > cymax)
891 else if (iy > cymin && iy < cymax && ybkge >= 1.0)
893 along_block(ix,iy,iz).y() +=
area;
894 betai_inv1(ix,iy,iz).y() += vol_block / (1.0 -
area +
floatSMALL);
896 else if (iy == cymin || iy == cymax)
900 const scalar double_f = (cymin == cymax ? 1.0 : 0.5);
902 add_blockage_c(area_block_s(ix,iy,iz).
y(), dirn_block(ix,iy,iz).
y(), area, double_f);
903 cell_area = (zgrid.
width(iz) * xgrid.
width(ix));
904 surf(ix,iy,iz) += double_f *
area * cell_area;
905 betai_inv1(ix,iy,iz).y() += vol_block / (1.0 -
area +
floatSMALL);
911 if (temp > 0.0) { grating_count(ix,iy,iz).y() += temp; }
915 area = xov * yov * zbkge;
916 if (iz < czmin || iz > czmax)
918 else if (iz > czmin && iz < czmax && zbkge >= 1.0)
920 along_block(ix,iy,iz).z() +=
area;
921 betai_inv1(ix,iy,iz).z() += vol_block / (1.0 -
area +
floatSMALL);
923 else if (iz == czmin || iz == czmax)
927 const scalar double_f = (czmin == czmax ? 1.0 : 0.5);
929 add_blockage_c(area_block_s(ix,iy,iz).z(), dirn_block(ix,iy,iz).z(), area, double_f);
930 cell_area = (xgrid.
width(ix) * ygrid.
width(iy));
931 surf(ix,iy,iz) += double_f *
area * cell_area;
932 betai_inv1(ix,iy,iz).z() += vol_block / (1.0 -
area +
floatSMALL);
938 if (temp > 0.0) { grating_count(ix,iy,iz).z() += temp; }
965 for (label ix = cxmin; ix <= cxmax; ix++)
968 scalar olap25 = 0.25 * xoverlap[ix];
971 ((
pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
973 obs_count(ix,cymin,czmin) += temp;
974 obs_count(ix,cymax,czmin) += temp;
975 obs_count(ix,cymin,czmax) += temp;
976 obs_count(ix,cymax,czmax) += temp;
978 sub_count(ix,cymin,czmin).x() += temp;
979 sub_count(ix,cymax,czmin).x() += temp;
980 sub_count(ix,cymin,czmax).x() += temp;
981 sub_count(ix,cymax,czmax).x() += temp;
986 drag_s(ix,cymin,czmin).yy() += zoverlap[czmin] * double_f.z() * olap25 * ybkge / ygrid.
width(cymin);
987 drag_s(ix,cymax,czmin).yy() += zoverlap[czmin] * double_f.z() * olap25 * ybkge / ygrid.
width(cymax);
988 drag_s(ix,cymin,czmax).yy() += zoverlap[czmax] * double_f.z() * olap25 * ybkge / ygrid.
width(cymin);
989 drag_s(ix,cymax,czmax).yy() += zoverlap[czmax] * double_f.z() * olap25 * ybkge / ygrid.
width(cymax);
991 drag_s(ix,cymin,czmin).zz() += yoverlap[cymin] * double_f.y() * olap25 * zbkge / zgrid.
width(czmin);
992 drag_s(ix,cymax,czmin).zz() += yoverlap[cymax] * double_f.y() * olap25 * zbkge / zgrid.
width(czmin);
993 drag_s(ix,cymin,czmax).zz() += yoverlap[cymin] * double_f.y() * olap25 * zbkge / zgrid.
width(czmax);
994 drag_s(ix,cymax,czmax).zz() += yoverlap[cymax] * double_f.y() * olap25 * zbkge / zgrid.
width(czmax);
999 for (label iy = cymin+1; iy < cymax; iy++)
1001 for (label iz = czmin+1; iz < czmax; iz++)
1004 drag_s(ix,iy,iz).xx() = xbkge / xgrid.
width(ix);
1010 for (label iy = cymin; iy <= cymax; iy++)
1012 scalar olap25 = 0.25 * yoverlap[iy];
1014 ((
pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
1016 obs_count(cxmin,iy,czmin) += temp;
1017 obs_count(cxmax,iy,czmin) += temp;
1018 obs_count(cxmin,iy,czmax) += temp;
1019 obs_count(cxmax,iy,czmax) += temp;
1021 sub_count(cxmin,iy,czmin).y() += temp;
1022 sub_count(cxmax,iy,czmin).y() += temp;
1023 sub_count(cxmin,iy,czmax).y() += temp;
1024 sub_count(cxmax,iy,czmax).y() += temp;
1028 if (iy > cymin && iy < cymax)
1030 drag_s(cxmin,iy,czmin).zz() += xoverlap[cxmin] * double_f.x() * olap25 * zbkge / zgrid.
width(czmin);
1031 drag_s(cxmax,iy,czmin).zz() += xoverlap[cxmin] * double_f.x() * olap25 * zbkge / zgrid.
width(czmin);
1032 drag_s(cxmin,iy,czmax).zz() += xoverlap[cxmax] * double_f.x() * olap25 * zbkge / zgrid.
width(czmax);
1033 drag_s(cxmax,iy,czmax).zz() += xoverlap[cxmax] * double_f.x() * olap25 * zbkge / zgrid.
width(czmax);
1035 drag_s(cxmin,iy,czmin).xx() += zoverlap[czmin] * double_f.z() * olap25 * xbkge / xgrid.
width(cxmin);
1036 drag_s(cxmax,iy,czmin).xx() += zoverlap[czmax] * double_f.z() * olap25 * xbkge / xgrid.
width(cxmax);
1037 drag_s(cxmin,iy,czmax).xx() += zoverlap[czmin] * double_f.z() * olap25 * xbkge / xgrid.
width(cxmin);
1038 drag_s(cxmax,iy,czmax).xx() += zoverlap[czmax] * double_f.z() * olap25 * xbkge / xgrid.
width(cxmax);
1043 for (label iz = czmin+1; iz < czmax; iz++)
1045 for (label ix = cxmin+1; ix < cxmax; ix++)
1048 drag_s(ix,iy,iz).yy() = ybkge / ygrid.
width(iy);
1054 for (label iz = czmin; iz <= czmax; iz++)
1056 scalar olap25 = 0.25 * zoverlap[iz];
1058 ((
pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
1060 obs_count(cxmin,cymin,iz) += temp;
1061 obs_count(cxmin,cymax,iz) += temp;
1062 obs_count(cxmax,cymin,iz) += temp;
1063 obs_count(cxmax,cymax,iz) += temp;
1065 sub_count(cxmin,cymin,iz).z() += temp;
1066 sub_count(cxmin,cymax,iz).z() += temp;
1067 sub_count(cxmax,cymin,iz).z() += temp;
1068 sub_count(cxmax,cymax,iz).z() += temp;
1072 if (iz > czmin && iz < czmax)
1074 drag_s(cxmin,cymin,iz).xx() += yoverlap[cymin] * double_f.y() * olap25 * xbkge / xgrid.
width(cxmin);
1075 drag_s(cxmax,cymin,iz).xx() += yoverlap[cymin] * double_f.y() * olap25 * xbkge / xgrid.
width(cxmax);
1076 drag_s(cxmin,cymax,iz).xx() += yoverlap[cymax] * double_f.y() * olap25 * xbkge / xgrid.
width(cxmin);
1077 drag_s(cxmax,cymax,iz).xx() += yoverlap[cymax] * double_f.y() * olap25 * xbkge / xgrid.
width(cxmax);
1079 drag_s(cxmin,cymin,iz).yy() += xoverlap[cxmin] * double_f.x() * olap25 * ybkge / ygrid.
width(cymin);
1080 drag_s(cxmax,cymin,iz).yy() += xoverlap[cxmax] * double_f.x() * olap25 * ybkge / ygrid.
width(cymin);
1081 drag_s(cxmin,cymax,iz).yy() += xoverlap[cxmin] * double_f.x() * olap25 * ybkge / ygrid.
width(cymax);
1082 drag_s(cxmax,cymax,iz).yy() += xoverlap[cxmax] * double_f.x() * olap25 * ybkge / ygrid.
width(cymax);
1088 for (label ix = cxmin+1; ix < cxmax; ix++)
1090 for (label iy = cymin+1; iy < cymax; iy++)
1093 drag_s(ix,iy,iz).zz() = zbkge / zgrid.
width(iz);
Various functions to operate on Lists.
Preparation of fields for PDRFoam.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
IjkField< vector > drag_r
Directional drag from round obstacles.
IjkField< vector > betai_inv1
IjkField< vector > along_block
Longitudinal area blockage from obstacles that extend all the way through the cell in a given directi...
IjkField< scalar > v_block
Volume blockage.
SquareMatrix< scalar > ac_lblock
void addBlockage(const PDRobstacle &obs, DynamicList< PDRpatchDef > &patches, const int volumeSign)
Add general (non-cylinder) blockage.
IjkField< vector > area_block_r
Summed area blockage (directional) from round obstacles.
Vector< List< scalar > > overlap_1d
SquareMatrix< scalar > b_lblock
IjkField< vector > sub_count
Number of obstacles parallel to specified direction.
IjkField< Vector< bool > > hole_in_face
Face field for (directional) hole in face.
SquareMatrix< symmTensor2D > c_drag
Cell-centred drag.
IjkField< vector > area_block_s
Summed area blockage (directional) from sharp obstacles.
SquareMatrix< scalar > abperim
SquareMatrix< scalar > bc_lblock
void addCylinder(const PDRobstacle &obs)
Add cylinder blockage.
const PDRblock & block() const
Reference to PDRblock.
SquareMatrix< scalar > aboverlap
IjkField< Vector< bool > > dirn_block
A total directional blockage in the cell.
IjkField< scalar > obs_count
Number of obstacles in cell.
IjkField< symmTensor > drag_s
Tensorial drag from sharp obstacles.
IjkField< vector > face_block
Face area blockage for face, summed from cell centre-plane to cell centre-plane.
SquareMatrix< scalar > a_lblock
IjkField< scalar > surf
Surface area in cell.
SquareMatrix< scalar > c_count
Grid locations in an axis direction.
scalar width(const label i) const
Cell size at element position.
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation....
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Obstacle definitions for PDR.
direction orient
The x/y/z orientation (0,1,2).
scalar x() const
Obstacle position accessors.
vector span
The obstacle dimensions (for boxes).
int typeId
The obstacle type-id.
scalar blockedFacePar
Faces with area blockage greater than this are blocked.
Bookkeeping for patch definitions.
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
label size() const noexcept
The number of entries in the list.
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
A class for handling words, derived from Foam::string.
static word validate(const std::string &s, const bool prefix=false)
Construct validated word (no invalid characters).
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
label find_if(const ListType &input, const UnaryPredicate &pred, const label start=0)
Find index of the first occurrence that satisfies the predicate.
Utilities for PDR (eg, for setFields).
void two_d_overlap(const UList< scalar > &a_olap, label amin, label amax, const UList< scalar > &b_olap, label bmin, label bmax, SquareMatrix< scalar > &ab_olap)
Combine two 1D overlaps.
void one_d_overlap(scalar xmin, scalar xmax, const PDRblock::location &grid, List< scalar > &olap, int *cmin, int *cmax, int *cfmin, int *cfmax)
Determine 1-D overlap locations for a geometric entity.
void circle_overlap(scalar ac, scalar bc, scalar dia, scalar theta, scalar wa, scalar wb, const PDRblock::location &agrid, label amin, label amax, const PDRblock::location &bgrid, label bmin, label bmax, SquareMatrix< scalar > &ab_olap, SquareMatrix< scalar > &ab_perim, SquareMatrix< scalar > &a_lblock, SquareMatrix< scalar > &ac_lblock, SquareMatrix< scalar > &c_count, SquareMatrix< symmTensor2D > &c_drag, SquareMatrix< scalar > &b_lblock, SquareMatrix< scalar > &bc_lblock)
Calculate the proportion of each (two-dimensional) grid cell overlapped by the circle or angled recta...
const wordList area
Standard area field types (scalar, vector, tensor, etc).
dimensionedScalar sin(const dimensionedScalar &ds)
bool equal(const T &a, const T &b)
Compare two values for equality.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Foam::PDRparams pars
Globals for program parameters (ugly hack).
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< scalar > scalarList
List of scalar.
dimensionedScalar cos(const dimensionedScalar &ds)
bool isNull(const T *ptr) noexcept
True if ptr is a pointer (of type T) to the nullObject.
constexpr char nl
The newline '\n' character (0x0a).