PDRarraysAnalyse.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) 2016 Shell Research Ltd.
9  Copyright (C) 2019 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 "PDRsetFields.H"
30 #include "PDRobstacle.H"
31 #include "PDRpatchDef.H"
32 #include "PDRutils.H"
33 #include "PDRutilsInternal.H"
34 #include "ListOps.H"
35 
36 #include <cstdlib>
37 #include <cstdio>
38 #include <cstring>
39 
40 #ifndef FULLDEBUG
41 #ifndef NDEBUG
42 #define NDEBUG
43 #endif
44 #endif
45 #include <cassert>
46 
47 using namespace Foam;
48 using namespace Foam::PDRutils;
49 
50 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
51 
52 // Cell blockage.
53 //
54 // Use simple sum, because this will eventually be divided by a notional
55 // number of rows to give a per-row blockage ratio
56 //
57 // b is the (fractional) area blockage. f is 1 or 0.5 to split between ends
58 // Thus if b isclose to 1, the obstacle is totally blocking the cell in this direction,
59 // and we could modify the behaviour if we wish.
60 inline static void add_blockage_c
61 (
62  scalar& a,
63  bool& blocked,
64  const scalar b,
65  const scalar f = 1.0
66 )
67 {
68  a += b * f;
69  if (b > pars.blockageNoCT)
70  {
71  blocked = true;
72  }
73 }
74 
75 
76 // Face blockage
77 //
78 // Adds more area blockage to existing amount by assuming partial overlap,
79 // i.e. multiplying porosities.
80 //
81 // Simple addition if the existing amount is negative, because negative
82 // blocks (summed first) should just cancel out part of positive blocks.
83 inline static void add_blockage_f
84 (
85  scalar& a,
86  const scalar b,
87  bool isHole
88 )
89 {
90  if (a > 0.0)
91  {
92  // Both positive
93  a = 1.0 - (1.0 - a) * (1.0 - b);
94  }
95  else if (b < pars.blockedFacePar || isHole)
96  {
97  // Add until it eventually becomes positive
98  a += b;
99  }
100  else
101  {
102  // If one obstacle blocks face, face is blocked, regardless of
103  // overlap calculations, unless an input negative obstacle makes a
104  // hole in it
105  a = b;
106  }
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  if (equal(obs.vbkge, 0))
115  {
116  return;
117  }
118 
119  if (isNull(block()))
120  {
122  << "No PDRblock set" << nl
123  << exit(FatalError);
124  }
125 
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();
130 
131  scalarList& xoverlap = overlap_1d.x();
132  scalarList& yoverlap = overlap_1d.y();
133  scalarList& zoverlap = overlap_1d.z();
134 
135  int cxmin, cxmax, cymin, cymax, czmin, czmax;
136  int cfxmin, cfxmax, cfymin, cfymax, cfzmin, cfzmax;
137 
138  scalar rad_a, rad_b;
139  vector area_block(Zero);
140 
141  if (obs.typeId == PDRobstacle::CYLINDER)
142  {
143  rad_a = rad_b = 0.5*obs.dia();
144  }
145  else
146  {
147  rad_a = 0.5*(obs.wa * ::cos(obs.theta()) + obs.wb * ::sin(obs.theta()));
148  rad_b = 0.5*(obs.wa * ::sin(obs.theta()) + obs.wb * ::cos(obs.theta()));
149  }
150 
151  switch (obs.orient)
152  {
153  case vector::Z:
154  {
155  // Determine the part of the grid (potentially) covered by this obstacle.
157  (
158  obs.x() - rad_a, obs.x() + rad_a,
159  pdrBlock.grid().x(),
160  xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
161  ); assert(cxmax >=0);
162 
164  (
165  obs.y() - rad_b, obs.y() + rad_b,
166  pdrBlock.grid().y(),
167  yoverlap, &cymin, &cymax, &cfymin, &cfymax
168  ); assert(cymax >=0);
169 
171  (
172  obs.z(), obs.z() + obs.len(),
173  pdrBlock.grid().z(),
174  zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
175  ); assert(czmax >=0);
176 
177  // The area of intersection with each 2D cell in an x-y plane.
178  // a corresponds to x, and b to y.
180  (
181  obs.x(), obs.y(), obs.dia(),
182  obs.theta(), obs.wa, obs.wb,
183  pdrBlock.grid().x(), cxmin, cfxmax,
184  pdrBlock.grid().y(), cymin, cfymax,
185  aboverlap, abperim, a_lblock, ac_lblock, c_count, c_drag,
186  b_lblock, bc_lblock
187  );
188 
189 
190  for (label ix = cxmin; ix <= cfxmax; ix++)
191  {
192  for (label iy = cymin; iy <= cfymax; iy++)
193  {
194  for (label iz = czmin; iz <= cfzmax; iz++)
195  {
196  const scalar vol_block = aboverlap(ix,iy) * zoverlap[iz];
197  v_block(ix,iy,iz) += vol_block;
198  surf(ix,iy,iz) += abperim(ix,iy) * zoverlap[iz] * zgrid.width(iz);
199 
200  // In the 2D case, the ends of the cylinder appear to
201  // be in the cell even when not, making surf and
202  // obs_size wrong. So leave out ends.
203 
204  if (!pars.two_d && (iz == czmin || iz == czmax))
205  {
206  // End cells
207  const scalar both_ends_fac = (czmin == czmax ? 2.0 : 1.0);
208 
209  surf(ix,iy,iz) += aboverlap(ix,iy)
210  * xgrid.width(ix) * ygrid.width(iy) * both_ends_fac;
211  }
212 
213  const scalar temp = c_count(ix,iy) * zoverlap[iz];
214 
215  obs_count(ix,iy,iz) += temp;
216  sub_count(ix,iy,iz).z() += temp;
217 
218  if (!pars.two_d && (iz == czmin || iz == czmax))
219  {
220  // End faces
221  const scalar both_ends_fac = (czmin == czmax ? 2.0 : 1.0);
222 
223  sub_count(ix,iy,iz).z() -= temp * both_ends_fac / 2.0;
224  }
225 
226  // Keep the blockage and drag of round obst separate
227  // from the sharp for the moment because different
228  // blockage ratio corrections will be needed later.
229  //
230  // Only the relevant diagonal elements of drag
231  // are stored; other are zero.
232 
233  area_block.x() = ac_lblock(ix,iy) * zoverlap[iz];
234  area_block.y() = bc_lblock(ix,iy) * zoverlap[iz];
235 
236  // Do not want blockage and drag across the end
237  // except at perimeter.
238  if (aboverlap(ix,iy) < pars.blockedFacePar)
239  {
240  if (obs.typeId == PDRobstacle::CYLINDER)
241  {
242  drag_r(ix,iy,iz).x() += c_drag(ix,iy).xx() * zoverlap[iz] / xgrid.width(ix) / ygrid.width(iy);
243  drag_r(ix,iy,iz).y() += c_drag(ix,iy).yy() * zoverlap[iz] / xgrid.width(ix) / ygrid.width(iy);
244 
245  add_blockage_c(area_block_r(ix,iy,iz).x(), dirn_block(ix,iy,iz).x(), area_block.x(), 1.0);
246  add_blockage_c(area_block_r(ix,iy,iz).y(), dirn_block(ix,iy,iz).y(), area_block.y(), 1.0);
247  }
248  else
249  {
250  drag_s(ix,iy,iz).xx() += c_drag(ix,iy).xx() * zoverlap[iz] / xgrid.width(ix) / ygrid.width(iy);
251  drag_s(ix,iy,iz).yy() += c_drag(ix,iy).yy() * zoverlap[iz] / xgrid.width(ix) / ygrid.width(iy);
252 
253  add_blockage_c(area_block_s(ix,iy,iz).x(), dirn_block(ix,iy,iz).x(), area_block.x(), 1.0);
254  add_blockage_c(area_block_s(ix,iy,iz).y(), dirn_block(ix,iy,iz).y(), area_block.y(), 1.0);
255  }
256  }
257  // Here we accumulate 1/betai - 1.
258  betai_inv1(ix,iy,iz).x() += vol_block / (1.0 - area_block.x() + floatSMALL);
259  betai_inv1(ix,iy,iz).y() += vol_block / (1.0 - area_block.y() + floatSMALL);
260  betai_inv1(ix,iy,iz).z() += vol_block / (1.0 - aboverlap(ix,iy) + floatSMALL);
261 
262  // The off-diagonal elements of drag are stored in "drag" for BOTH round and sharp ostacles
263  drag_s(ix,iy,iz).xy() += c_drag(ix,iy).xy() * zoverlap[iz] / xgrid.width(ix) / ygrid.width(iy);
264 
265  add_blockage_f(face_block(ix,iy,iz).x(), a_lblock(ix,iy) * zoverlap[iz], hole_in_face(ix,iy,iz).x());
266  add_blockage_f(face_block(ix,iy,iz).y(), b_lblock(ix,iy) * zoverlap[iz], hole_in_face(ix,iy,iz).y());
267  }
268 
269  // Face blockage in the z direction
270  if (czmin == czmax)
271  {
272  // Does not span cell. Block nearest face.
273  add_blockage_f(face_block(ix,iy,cfzmax).z(), aboverlap(ix,iy), hole_in_face(ix,iy,cfzmax).z());
274  }
275  else
276  {
277  // In at least two cells.
278  // Block first and last faces overlapped
279  add_blockage_f(face_block(ix,iy,czmin+1).z(), aboverlap(ix,iy), hole_in_face(ix,iy,czmin+1).z());
280  if (czmax > czmin+1)
281  {
282  add_blockage_f(face_block(ix,iy,czmax).z(), aboverlap(ix,iy), hole_in_face(ix,iy,czmax).z());
283  }
284  }
285 
286  // z_block is used to work out the blockage ratio for each
287  // "row" of sub-grid obstacles so this cylinder should
288  // not have any eeffct in cells that it completely spans.
289  // Hence statement commented out below and replaced by
290  // two lines after this loop. That longitudinal clockage
291  // goes into new array z_aling_block
292 
293  for (label iz = czmin+1; iz < czmax; ++iz)
294  {
295  // Internal only
296  along_block(ix,iy,iz).z() += aboverlap(ix,iy);
297  }
298 
299  // Longitudinal drag only applies at ends of cylinder.
300  // If cylinder spans more than 1 cell, apply half at each
301  // end.
302 
303  drag_s(ix,iy,czmin).zz() += aboverlap(ix,iy) * xgrid.width(ix) * ygrid.width(iy) / 2.0;
304  drag_s(ix,iy,czmax).zz() += aboverlap(ix,iy) * xgrid.width(ix) * ygrid.width(iy) / 2.0;
305 
306  add_blockage_c(area_block_s(ix,iy,czmin).z(), dirn_block(ix,iy,czmin).z(), aboverlap(ix,iy), 0.5);
307  add_blockage_c(area_block_s(ix,iy,czmax).z(), dirn_block(ix,iy,czmax).z(), aboverlap(ix,iy), 0.5);
308  }
309  }
310 
311  break;
312  }
313 
314  case vector::X: // orientation
315  {
316  // x-direction cylinder. a,b are y,z.
318  (
319  obs.x(), obs.x() + obs.len(),
320  pdrBlock.grid().x(),
321  xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
322  ); assert(cxmax >=0);
323 
325  (
326  obs.y() - rad_a, obs.y() + rad_a,
327  pdrBlock.grid().y(),
328  yoverlap, &cymin, &cymax, &cfymin, &cfymax
329  ); assert(cymax >=0);
330 
332  (
333  obs.z() - rad_b, obs.z() + rad_b,
334  pdrBlock.grid().z(),
335  zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
336  ); assert(czmax >=0);
337 
339  (
340  obs.y(), obs.z(), obs.dia(),
341  obs.theta(), obs.wa, obs.wb,
342  pdrBlock.grid().y(), cymin, cfymax,
343  pdrBlock.grid().z(), czmin, cfzmax,
344  aboverlap, abperim, a_lblock, ac_lblock, c_count, c_drag,
345  b_lblock, bc_lblock
346  );
347 
348 
349  for (label iy = cymin; iy <= cfymax; iy++)
350  {
351  for (label iz = czmin; iz <= cfzmax; iz++)
352  {
353  for (label ix = cxmin; ix <= cxmax; ix++)
354  {
355  const scalar vol_block = aboverlap(iy,iz) * xoverlap[ix];
356  v_block(ix,iy,iz) += vol_block;
357  surf(ix,iy,iz) += abperim(iy,iz) * xoverlap[ix] * xgrid.width(ix);
358 
359  if (ix == cxmin || ix == cxmax)
360  {
361  // End cells
362  const scalar both_ends_fac = (cxmin == cxmax ? 2.0 : 1.0);
363 
364  surf(ix,iy,iz) += aboverlap(iy,iz)
365  * ygrid.width(iy) * zgrid.width(iz) * both_ends_fac;
366  }
367 
368  const scalar temp = c_count(iy,iz) * xoverlap[ix];
369  obs_count(ix,iy,iz) += temp;
370  sub_count(ix,iy,iz).x() += temp;
371 
372  if (ix == cfxmin || ix == cfxmax)
373  {
374  // End faces
375  const scalar both_ends_fac = (cfxmin == cfxmax ? 2.0 : 1.0);
376 
377  sub_count(ix,iy,iz).x() -= temp * both_ends_fac / 2.0;
378  }
379 
380  area_block.y() = ac_lblock(iy,iz) * xoverlap[ix];
381  area_block.z() = bc_lblock(iy,iz) * xoverlap[ix];
382 
383  // Do not want blockage and drag across the end
384  // except at perimeter.
385  if (aboverlap(iy,iz) < pars.blockedFacePar)
386  {
387  if (obs.typeId == PDRobstacle::CYLINDER)
388  {
389  drag_r(ix,iy,iz).y() += c_drag(iy,iz).xx() * xoverlap[ix] / ygrid.width(iy) / zgrid.width(iz);
390  drag_r(ix,iy,iz).z() += c_drag(iy,iz).yy() * xoverlap[ix] / ygrid.width(iy) / zgrid.width(iz);
391 
392  add_blockage_c(area_block_r(ix,iy,iz).y(), dirn_block(ix,iy,iz).y(), area_block.y(), 1.0);
393  add_blockage_c(area_block_r(ix,iy,iz).z(), dirn_block(ix,iy,iz).z(), area_block.z(), 1.0);
394  }
395  else
396  {
397  drag_s(ix,iy,iz).yy() += c_drag(iy,iz).xx() * xoverlap[ix] / ygrid.width(iy) / zgrid.width(iz);
398  drag_s(ix,iy,iz).zz() += c_drag(iy,iz).yy() * xoverlap[ix] / ygrid.width(iy) / zgrid.width(iz);
399 
400  add_blockage_c(area_block_s(ix,iy,iz).y(), dirn_block(ix,iy,iz).y(), area_block.y(), 1.0);
401  add_blockage_c(area_block_s(ix,iy,iz).z(), dirn_block(ix,iy,iz).z(), area_block.z(), 1.0);
402  }
403  }
404  betai_inv1(ix,iy,iz).y() += vol_block / (1.0 - area_block.y() + floatSMALL);
405  betai_inv1(ix,iy,iz).z() += vol_block / (1.0 - area_block.z() + floatSMALL);
406  betai_inv1(ix,iy,iz).x() += vol_block / (1.0 - aboverlap(iy,iz) + floatSMALL);
407 
408  // The off-diagonal elements of drag are stored in "drag" for BOTH round and sharp ostacles
409  drag_s(ix,iy,iz).yz() += c_drag(iy,iz).xy() * xoverlap[ix] / ygrid.width(iy) / zgrid.width(iz);
410 
411  add_blockage_f(face_block(ix,iy,iz).y(), a_lblock(iy,iz) * xoverlap[ix], hole_in_face(ix,iy,iz).y());
412  add_blockage_f(face_block(ix,iy,iz).z(), b_lblock(iy,iz) * xoverlap[ix], hole_in_face(ix,iy,iz).z());
413  }
414  if (cxmin == cxmax)
415  {
416  // Does not span cell. Block nearest face.
417  add_blockage_f(face_block(cfxmax,iy,iz).x(), aboverlap(iy,iz), hole_in_face(cfxmax,iy,iz).x());
418  }
419  else
420  {
421  // In at least two cells.
422  // Block first and last faces overlapped
423  add_blockage_f(face_block(cxmin+1,iy,iz).x(), aboverlap(iy,iz), hole_in_face(cxmin+1,iy,iz).x());
424  if (cxmax > cxmin+1)
425  {
426  add_blockage_f(face_block(cxmax,iy,iz).x(), aboverlap(iy,iz), hole_in_face(cxmax,iy,iz).x());
427  }
428  }
429 
430  for (label ix = cxmin+1; ix < cxmax; ++ix)
431  {
432  // Internal only
433  along_block(ix,iy,iz).x() += aboverlap(iy,iz);
434  }
435  drag_s(cxmin,iy,iz).xx() += aboverlap(iy,iz) * ygrid.width(iy) * zgrid.width(iz) / 2.0;
436  drag_s(cxmax,iy,iz).xx() += aboverlap(iy,iz) * ygrid.width(iy) * zgrid.width(iz) / 2.0;
437 
438  add_blockage_c(area_block_s(cxmin,iy,iz).x(), dirn_block(cxmin,iy,iz).x(), aboverlap(iy,iz), 0.5);
439  add_blockage_c(area_block_s(cxmax,iy,iz).x(), dirn_block(cxmax,iy,iz).x(), aboverlap(iy,iz), 0.5);
440  }
441  }
442 
443  break;
444  }
445 
446  case vector::Y: // orientation
447  {
448  // y-direction cylinder. a,b are z,x.
450  (
451  obs.x() - rad_b, obs.x() + rad_b,
452  pdrBlock.grid().x(),
453  xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
454  ); assert(cxmax >=0);
455 
457  (
458  obs.y(), obs.y() + obs.len(),
459  pdrBlock.grid().y(),
460  yoverlap, &cymin, &cymax, &cfymin, &cfymax
461  ); assert(cymax >=0);
462 
464  (
465  obs.z() - rad_a, obs.z() + rad_a,
466  pdrBlock.grid().z(),
467  zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
468  ); assert(czmax >=0);
469 
471  (
472  obs.z(), obs.x(), obs.dia(),
473  obs.theta(), obs.wa, obs.wb,
474  pdrBlock.grid().z(), czmin, cfzmax,
475  pdrBlock.grid().x(), cxmin, cfxmax,
476  aboverlap, abperim, a_lblock, ac_lblock, c_count, c_drag,
477  b_lblock, bc_lblock
478  );
479 
480 
481  for (label iz = czmin; iz <= cfzmax; iz++)
482  {
483  for (label ix = cxmin; ix <= cfxmax; ix++)
484  {
485  for (label iy = cymin; iy <= cymax; iy++)
486  {
487  const scalar vol_block = aboverlap(iz,ix) * yoverlap[iy];
488  v_block(ix,iy,iz) += vol_block;
489  surf(ix,iy,iz) += abperim(iz,ix) * yoverlap[iy] * ygrid.width(iy);
490 
491  if (iy == cymin || iy == cymax)
492  {
493  // End cells
494  const scalar both_ends_fac = (cymin == cymax ? 2.0 : 1.0);
495 
496  surf(ix,iy,iz) += aboverlap(iz,ix)
497  * zgrid.width(iz) * xgrid.width(ix) * both_ends_fac;
498  }
499 
500  const scalar temp = c_count(iz,ix) * yoverlap[iy];
501 
502  obs_count(ix,iy,iz) += temp;
503  sub_count(ix,iy,iz).y() += temp;
504 
505  if (iy == cfymin || iy == cfymax)
506  {
507  // End faces
508  const scalar both_ends_fac = (cfymin == cfymax ? 2.0 : 1.0);
509 
510  sub_count(ix,iy,iz).y() -= temp * both_ends_fac / 2.0;
511  }
512 
513  area_block.z() = ac_lblock(iz,ix) * yoverlap[iy];
514  area_block.x() = bc_lblock(iz,ix) * yoverlap[iy];
515 
516  // Do not want blockage and drag across the end
517  // except at perimeter.
518  if (aboverlap(iz,ix) < pars.blockedFacePar)
519  {
520  if (obs.typeId == PDRobstacle::CYLINDER)
521  {
522  drag_r(ix,iy,iz).z() += c_drag(iz,ix).xx() * yoverlap[iy] / zgrid.width(iz) / xgrid.width(ix);
523  drag_r(ix,iy,iz).x() += c_drag(iz,ix).yy() * yoverlap[iy] / zgrid.width(iz) / xgrid.width(ix);
524 
525  add_blockage_c(area_block_r(ix,iy,iz).z(), dirn_block(ix,iy,iz).z(), area_block.z(), 1.0);
526  add_blockage_c(area_block_r(ix,iy,iz).x(), dirn_block(ix,iy,iz).x(), area_block.x(), 1.0);
527  }
528  else
529  {
530  drag_s(ix,iy,iz).zz() += c_drag(iz,ix).xx() * yoverlap[iy] / zgrid.width(iz) / xgrid.width(ix);
531  drag_s(ix,iy,iz).xx() += c_drag(iz,ix).yy() * yoverlap[iy] / zgrid.width(iz) / xgrid.width(ix);
532 
533  add_blockage_c(area_block_s(ix,iy,iz).z(), dirn_block(ix,iy,iz).z(), area_block.z(), 1.0);
534  add_blockage_c(area_block_s(ix,iy,iz).x(), dirn_block(ix,iy,iz).x(), area_block.x(), 1.0);
535  }
536  }
537  betai_inv1(ix,iy,iz).z() += vol_block / (1.0 - area_block.z() + floatSMALL);
538  betai_inv1(ix,iy,iz).x() += vol_block / (1.0 - area_block.x() + floatSMALL);
539  betai_inv1(ix,iy,iz).y() += vol_block / (1.0 - aboverlap(iz,ix) + floatSMALL);
540 
541  // The off-diagonal elements of drag are stored in "drag" for BOTH round and sharp obstacles
542  drag_s(ix,iy,iz).xz() += c_drag(iz,ix).xy() * yoverlap[iy] / zgrid.width(iz) / xgrid.width(ix);
543 
544  add_blockage_f(face_block(ix,iy,iz).z(), a_lblock(iz,ix) * yoverlap[iy], hole_in_face(ix,iy,iz).z());
545  add_blockage_f(face_block(ix,iy,iz).x(), b_lblock(iz,ix) * yoverlap[iy], hole_in_face(ix,iy,iz).x());
546  }
547  if (cymin == cymax)
548  {
549  // Does not span cell. Block nearest face.
550  add_blockage_f(face_block(ix,cfymax,iz).y(), aboverlap(iz,ix), hole_in_face(ix,cfymax,iz).y());
551  }
552  else
553  {
554  // In at least two cells.
555  // Block first and last faces overlapped
556  add_blockage_f(face_block(ix,cymin+1,iz).y(), aboverlap(iz,ix), hole_in_face(ix,cymin+1,iz).y());
557  if (cymax > cymin+1)
558  {
559  add_blockage_f(face_block(ix,cymax,iz).y(), aboverlap(iz,ix), hole_in_face(ix,cymax,iz).y());
560  }
561  }
562 
563  for (label iy = cymin+1; iy < cymax; ++iy)
564  {
565  // Internal only
566  along_block(ix,iy,iz).y() += aboverlap(iz,ix);
567  }
568 
569  drag_s(ix,cymin,iz).yy() += aboverlap(iz,ix) * zgrid.width(iz) * xgrid.width(ix) / 2.0;
570  drag_s(ix,cymax,iz).yy() += aboverlap(iz,ix) * zgrid.width(iz) * xgrid.width(ix) / 2.0;
571 
572  add_blockage_c(area_block_s(ix,cymin,iz).y(), dirn_block(ix,cymin,iz).y(), aboverlap(iz,ix), 0.5);
573  add_blockage_c(area_block_s(ix,cymax,iz).y(), dirn_block(ix,cymax,iz).y(), aboverlap(iz,ix), 0.5);
574  }
575  }
576 
577  break;
578  }
579 
580  default: // orientation
581  {
582  Info<< "Unexpected orientation " << int(obs.orient) << nl;
583  break;
584  }
585  }
586 }
587 
588 
590 (
591  const PDRobstacle& obs,
593  const int volumeSign
594 )
595 {
596  // The volumeSign indicates whether this pass is for negative or positive
597  // obstacles. Both if 0.
598  if
599  (
600  equal(obs.vbkge, 0)
601  || (volumeSign < 0 && obs.vbkge >= 0)
602  || (volumeSign > 0 && obs.vbkge < 0)
603  )
604  {
605  return;
606  }
607 
608  if (isNull(block()))
609  {
611  << "No PDRblock set" << nl
612  << exit(FatalError);
613  }
614 
615  const PDRblock& pdrBlock = block();
616  const PDRblock::location& xgrid = pdrBlock.grid().x();
617  const PDRblock::location& ygrid = pdrBlock.grid().y();
618  const PDRblock::location& zgrid = pdrBlock.grid().z();
619 
620  scalarList& xoverlap = overlap_1d.x();
621  scalarList& yoverlap = overlap_1d.y();
622  scalarList& zoverlap = overlap_1d.z();
623 
624 
625  // 0 will be used later for faces found to be blocked.
626  // 2 is used for wall-function faces.
627  label patchNum = PDRpatchDef::LAST_PREDEFINED;
628 
629  // Only the part of the panel that covers full cell faces will be used
630  // so later should keep the panel in the list plus a hole (-ve obstacle) for
631  // part that will become blowoff b.c.
632 
633  int indir = 0;
634 
635  // Panel or patch
636  const bool isPatch =
637  (
639  || (obs.typeId == PDRobstacle::RECT_PATCH)
640  );
641 
642  if (isPatch)
643  {
644  const auto& identifier = obs.identifier;
645 
646  const auto spc = identifier.find_first_of(" \t\n\v\f\r");
647 
648  const word patchName = word::validate(identifier.substr(0, spc));
649 
650  patchNum = ListOps::find
651  (
652  patches,
653  [=](const PDRpatchDef& p){ return patchName == p.patchName; },
654  1 // skip 0 (blocked face)
655  );
656 
657  if (patchNum < 1)
658  {
659  // The patch name was not already in the list
660  patchNum = patches.size();
661 
662  patches.append(PDRpatchDef(patchName));
663  }
664 
665 
666  PDRpatchDef& p = patches[patchNum];
667 
668  if (obs.typeId == PDRobstacle::RECT_PATCH)
669  {
670  indir = obs.inlet_dirn;
671  p.patchType = 0;
672  }
673  else
674  {
675  p.patchType = obs.blowoff_type;
676  p.blowoffPress = obs.blowoff_press;
677  p.blowoffTime = obs.blowoff_time;
678  if (obs.span.x() < 0.01)
679  {
680  indir = 1;
681  }
682  else if (obs.span.y() < 0.01)
683  {
684  indir = 2;
685  }
686  else if (obs.span.z() < 0.01)
687  {
688  indir = 3;
689  }
690  else
691  {
693  << "Blowoff panel should have zero thickness" << nl
694  << exit(FatalError);
695  }
696  }
697  }
698 
699  int cxmin, cxmax, cfxmin, cfxmax;
701  (
702  obs.x(), obs.x() + obs.span.x(),
703  pdrBlock.grid().x(),
704  xoverlap, &cxmin, &cxmax, &cfxmin, &cfxmax
705  ); assert(cxmax >=0);
706 
707  int cymin, cymax, cfymin, cfymax;
709  (
710  obs.y(), obs.y() + obs.span.y(),
711  pdrBlock.grid().y(),
712  yoverlap, &cymin, &cymax, &cfymin, &cfymax
713  ); assert(cymax >=0);
714 
715  int czmin, czmax, cfzmin, cfzmax;
717  (
718  obs.z(), obs.z() + obs.span.z(),
719  pdrBlock.grid().z(),
720  zoverlap, &czmin, &czmax, &cfzmin, &cfzmax
721  ); assert(czmax >=0);
722 
723  two_d_overlap(xoverlap, cxmin, cxmax, yoverlap, cymin, cymax, aboverlap);
724 
725 
726  const scalar vbkge = obs.vbkge;
727  const scalar xbkge = obs.xbkge;
728  const scalar ybkge = obs.ybkge;
729  const scalar zbkge = obs.zbkge;
730 
731  // Compensate for double-counting of drag if two edges in same cell
732  const vector double_f
733  (
734  ((cxmin == cxmax) ? 0.5 : 1.0),
735  ((cymin == cymax) ? 0.5 : 1.0),
736  ((czmin == czmax) ? 0.5 : 1.0)
737  );
738 
739  for (label ix = cxmin; ix <= cfxmax; ix++)
740  {
741  const scalar xov = xoverlap[ix];
742 
743  scalar area, cell_area, temp;
744 
745  for (label iy = cymin; iy <= cfymax; iy++)
746  {
747  const scalar yov = yoverlap[iy];
748 
749  for (label iz = czmin; iz <= cfzmax; iz++)
750  {
751  const scalar zov = zoverlap[iz];
752 
753  if
754  (
755  isPatch
756  &&
757  (
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)
764  )
765  )
766  {
767  /* Type RECT_PATCH (16) exists to set all faces it covers to be in a particular patch
768  usually an inlet or outlet.
769  ?? If the face not on a cell boundary, this will move it to the lower-cordinate
770  face of the relevant cell. It should be at the face of teh volume blocked by
771  the obstacle. But, if that is not at a cell boundary, the obstacle will be
772  putting blockage in front of teh vent, so we should be checking that it is
773  at a cell boundary. */
774 
775  switch (indir) // Face orientation
776  {
777  // X
778  case -1:
779  case 1:
780  if (yov * zov > pars.blockedFacePar)
781  {
782  face_patch(ix,iy,iz).x() = patchNum;
783  }
784  break;
785 
786  // Y
787  case -2:
788  case 2:
789  if (zov * xov > pars.blockedFacePar)
790  {
791  face_patch(ix,iy,iz).y() = patchNum;
792  }
793  break;
794 
795  // Z
796  case -3:
797  case 3:
798  if (xov * yov > pars.blockedFacePar)
799  {
800  face_patch(ix,iy,iz).z() = patchNum;
801  }
802  break;
803  }
804  } // End of code for patch
805 
806 
807  const scalar vol_block = aboverlap(ix,iy) * zov * vbkge;
808  v_block(ix,iy,iz) += vol_block;
809 
810  // These are the face blockages
811  if ((ix > cxmin && ix <= cxmax) || (cxmin == cxmax && ix == cfxmax))
812  {
813  temp = yov * zov * xbkge;
814 
815  // Has -ve volumeSign only when processing user-defined
816  // -ve blocks
817  if (volumeSign < 0)
818  {
819  hole_in_face(ix,iy,iz).x() = true;
820  }
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)
823  {
824  // Put faces of block in wall patch
825  face_patch(ix,iy,iz).x() = PDRpatchDef::WALL_PATCH;
826  }
827  }
828  if ((iy > cymin && iy <= cymax) || (cymin == cymax && iy == cfymax))
829  {
830  temp = zov * xov * ybkge;
831  if (volumeSign < 0)
832  {
833  hole_in_face(ix,iy,iz).y() = true;
834  }
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)
837  {
838  face_patch(ix,iy,iz).y() = PDRpatchDef::WALL_PATCH;
839  }
840  }
841  if ((iz > czmin && iz <= czmax) || (czmin == czmax && iz == cfzmax))
842  {
843  temp = xov * yov * zbkge;
844  if (volumeSign < 0)
845  {
846  hole_in_face(ix,iy,iz).z() = true;
847  }
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)
850  {
851  face_patch(ix,iy,iz).z() = PDRpatchDef::WALL_PATCH;
852  }
853  }
854 
855  // These are the interior blockages
856  /* As for cylinders, blockage that extends longitudinally all the way through the cell
857  should not be in x_block etc., but it does go into new arrays along_block etc. */
858  area = yov * zov * xbkge; // Note this is fraction of cell c-s area
859  if (ix < cxmin || ix > cxmax)
860  {}
861  else if (ix > cxmin && ix < cxmax && xbkge >= 1.0)
862  {
863  along_block(ix,iy,iz).x() += area;
864  betai_inv1(ix,iy,iz).x() += vol_block / (1.0 - area + floatSMALL);
865  }
866  else if (ix == cxmin || ix == cxmax)
867  {
868  // If front and back of the obstacle are not in the
869  // same cell, put half in each
870  const scalar double_f = (cxmin == cxmax ? 1.0 : 0.5);
871 
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);
876 
877  // To get Lobs right for a grating, allow for the fact that it is series of bars by having additional count
878  if (obs.typeId == PDRobstacle::GRATING && obs.orient == vector::X)
879  {
880  // * cell_area to make dimensional, then / std::sqrt(cell_area) to get width
881  temp = area * std::sqrt(cell_area) / obs.slat_width - 1.0;
882  if (temp > 0.0) { grating_count(ix,iy,iz).x() += temp; }
883  }
884  }
885 
886  /* As for cylinders, blockage that extends longitudinally all the way through the cell
887  should not be in x_block etc., but it does go into new arrays along_block etc. */
888  area = zov * xov * ybkge;
889  if (iy < cymin || iy > cymax)
890  {}
891  else if (iy > cymin && iy < cymax && ybkge >= 1.0)
892  {
893  along_block(ix,iy,iz).y() += area;
894  betai_inv1(ix,iy,iz).y() += vol_block / (1.0 - area + floatSMALL);
895  }
896  else if (iy == cymin || iy == cymax)
897  {
898  // If front and back of the obstacle are not in the
899  // same cell, put half in each
900  const scalar double_f = (cymin == cymax ? 1.0 : 0.5);
901 
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);
906 
907  if (obs.typeId == PDRobstacle::GRATING && obs.orient == vector::Y)
908  {
909  // * cell_area to make dimensional, then / std::sqrt(cell_area) to get width
910  temp = area * std::sqrt(cell_area) / obs.slat_width - 1.0;
911  if (temp > 0.0) { grating_count(ix,iy,iz).y() += temp; }
912  }
913  }
914 
915  area = xov * yov * zbkge;
916  if (iz < czmin || iz > czmax)
917  {}
918  else if (iz > czmin && iz < czmax && zbkge >= 1.0)
919  {
920  along_block(ix,iy,iz).z() += area;
921  betai_inv1(ix,iy,iz).z() += vol_block / (1.0 - area + floatSMALL);
922  }
923  else if (iz == czmin || iz == czmax)
924  {
925  // If front and back of the obstacle are not in the
926  // same cell, put half in each
927  const scalar double_f = (czmin == czmax ? 1.0 : 0.5);
928 
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);
933 
934  if (obs.typeId == PDRobstacle::GRATING && obs.orient == vector::Z)
935  {
936  // * cell_area to make dimensional, then / std::sqrt(cell_area) to get width
937  temp = area * std::sqrt(cell_area) / obs.slat_width - 1.0;
938  if (temp > 0.0) { grating_count(ix,iy,iz).z() += temp; }
939  }
940  }
941  }
942  }
943  }
944 
945  if (obs.typeId == PDRobstacle::RECT_PATCH)
946  {
947  // Was only needed to set face_patch values
948  return;
949  }
950 
951 
952  /* A narrow obstacle completely crossing the cell adds 1 to the count for the transverse directions
953  If all four edges are not in the relevant cell, we take 1/4 for each edge that is in.
954  If it does not totally cross the cell, then the count is reduced proportionately
955  If it is porous, then the count is reduced proportionately.
956  ?? Should it be? At least this is consistent with effect going smoothly to zero as porosity approaches 1 ??
957 
958  Note that more than 1 can be added for one obstacle, if sides and ends are in the cell.
959 
960  We do all the x-aligned edges first, both for y-blockage and z-blockage. */
961 
962  /* Intersection obstacles can have more edges than the intersecting obstacles that
963  generated them, so not good to incorporate these in N and drag. */
964 
965  for (label ix = cxmin; ix <= cxmax; ix++)
966  {
967  // Factor of 0.25 because 4 edges to be done
968  scalar olap25 = 0.25 * xoverlap[ix];
969 
970  const scalar temp =
971  ((pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
972 
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;
977 
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;
982 
983  // The 0.25 becomes 0.5 to allow for front/back faces in drag direction
984  olap25 *= 2.0;
985 
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);
990 
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);
995 
996  // Porous obstacles do not only have drag at edges
997  if (xbkge < 1.0)
998  {
999  for (label iy = cymin+1; iy < cymax; iy++)
1000  {
1001  for (label iz = czmin+1; iz < czmax; iz++)
1002  {
1003  // Internal only
1004  drag_s(ix,iy,iz).xx() = xbkge / xgrid.width(ix);
1005  }
1006  }
1007  }
1008  }
1009 
1010  for (label iy = cymin; iy <= cymax; iy++)
1011  {
1012  scalar olap25 = 0.25 * yoverlap[iy];
1013  const scalar temp =
1014  ((pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
1015 
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;
1020 
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;
1025 
1026  olap25 *= 2.0;
1027 
1028  if (iy > cymin && iy < cymax) // Avoid re-doing corners already done above
1029  {
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);
1034  }
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);
1039 
1040  // Porous obstacles do not only have drag at edges
1041  if (ybkge < 1.0)
1042  {
1043  for (label iz = czmin+1; iz < czmax; iz++)
1044  {
1045  for (label ix = cxmin+1; ix < cxmax; ix++)
1046  {
1047  // Internal only
1048  drag_s(ix,iy,iz).yy() = ybkge / ygrid.width(iy);
1049  }
1050  }
1051  }
1052  }
1053 
1054  for (label iz = czmin; iz <= czmax; iz++)
1055  {
1056  scalar olap25 = 0.25 * zoverlap[iz];
1057  const scalar temp =
1058  ((pars.noIntersectN && vbkge < 0.334) ? 0 : (olap25 * vbkge));
1059 
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;
1064 
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;
1069 
1070  olap25 *= 2.0;
1071 
1072  if (iz > czmin && iz < czmax) // Avoid re-doing corners already done above
1073  {
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);
1078 
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);
1083  }
1084 
1085  // Porous obstacles do not only have drag at edges
1086  if (zbkge < 1.0)
1087  {
1088  for (label ix = cxmin+1; ix < cxmax; ix++)
1089  {
1090  for (label iy = cymin+1; iy < cymax; iy++)
1091  {
1092  // Internal only
1093  drag_s(ix,iy,iz).zz() = zbkge / zgrid.width(iz);
1094  }
1095  }
1096  }
1097  }
1098 }
1099 
1100 
1101 // ************************************************************************* //
static word validate(const std::string &s, const bool prefix=false)
Construct validated word (no invalid characters).
Definition: word.C:39
scalar blockedFacePar
Faces with area blockage greater than this are blocked.
Definition: PDRparams.H:148
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.
label find(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as ListOps::find_if.
Definition: ListOps.H:790
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.
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...
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition: label.H:164
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
scalar blockageNoCT
If a single obstacle blocks a cell by more than this, then no CT in that direction.
Definition: PDRparams.H:159
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Utilities for PDR (eg, for setFields)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Various functions to operate on Lists.
scalar y
bool noIntersectN
Definition: PDRparams.H:81
scalar z() const
Definition: PDRobstacle.H:257
dimensionedScalar cos(const dimensionedScalar &ds)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
scalar y() const
Definition: PDRobstacle.H:256
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
const wordList area
Standard area field types (scalar, vector, tensor, etc)
A class for handling words, derived from Foam::string.
Definition: word.H:63
scalar dia() const
Definition: PDRobstacle.H:144
scalar len() const
Definition: PDRobstacle.H:146
Grid locations in an axis direction.
Definition: PDRblock.H:181
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation. Some of the input is similar to blockMeshDict, but since this specialization is for a single-block that is aligned with the x-y-z directions, it provides a different means of specifying the mesh.
Definition: PDRblock.H:149
Preparation of fields for PDRFoam.
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...
#define floatSMALL
Definition: PDRsetFields.H:48
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:54
dimensionedScalar sin(const dimensionedScalar &ds)
bool isNull(const T *ptr)
True if ptr is a pointer (of type T) to the nullObject.
Definition: nullObject.H:227
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
labelList f(nPoints)
scalar width(const label i) const
Cell size at element position.
Definition: PDRblockI.H:84
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition: PtrList.H:323
PtrList< volScalarField > & Y
vector span
The obstacle dimensions (for boxes)
Definition: PDRobstacle.H:140
const polyBoundaryMesh & patches
Foam::PDRparams pars
Globals for program parameters (ugly hack)
messageStream Info
Information stream (stdout output on master, null elsewhere)
int typeId
The obstacle type-id.
Definition: PDRobstacle.H:118
void addCylinder(const PDRobstacle &obs)
Add cylinder blockage.
volScalarField & p
direction orient
The x/y/z orientation (0,1,2)
Definition: PDRobstacle.H:123
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:49
void addBlockage(const PDRobstacle &obs, DynamicList< PDRpatchDef > &patches, const int volumeSign)
Add general (non-cylinder) blockage.
scalar x() const
Obstacle position accessors.
Definition: PDRobstacle.H:255
Obstacle definitions for PDR.
Definition: PDRobstacle.H:70
Namespace for OpenFOAM.
scalar theta() const
Definition: PDRobstacle.H:145
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Definition: PDRblock.H:711
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133