Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
embed-tree.h
Go to the documentation of this file.
1/** @file embed-tree.h
2 */
3/**
4# Embedded boundaries on adaptive trees
5
6This file defines the restriction/prolongation functions which are
7necessary to implement [embedded boundaries](embed.h) on adaptive
8meshes.
9
10## Volume fraction field *cs*
11
12For the embedded fraction field *cs*, the function below is modelled
13closely on the volume fraction refinement function
14[fraction_refine()](fractions.h#fraction_refine). */
15
17{
18 double cc = cs[];
19
20 /**
21 If the cell is empty or full, simple injection from the coarse cell
22 value is used. */
23
24 if (cc <= 0. || cc >= 1.) {
25 for (int _c = 0; _c < 4; _c++) /* foreach_child */
26 cs[] = cc;
27 }
28 else {
29
30 /**
31 If the cell contains the embedded boundary, we reconstruct the
32 boundary using VOF linear reconstruction and a normal estimated
33 from the surface fractions. */
34
36 double alpha = plane_alpha (cc, n);
37
38 for (int _c = 0; _c < 4; _c++) /* foreach_child */ {
39 static const coord a = {0.,0.,0.}, b = {.5,.5,.5};
40 coord nc;
41 for (int _d = 0; _d < dimension; _d++)
42 nc.x = child.x*n.x;
44 }
45 }
46}
47
48/**
49## Surface fractions field *fs*
50
51The embedded surface fractions *fs* are reconstructed using this
52function. */
53
54for (int _d = 0; _d < dimension; _d++)
56{
57 vector fs = s.v;
58
59 /**
60 If the cell is empty or full, simple injection from the coarse cell
61 value is used. */
62
63 if (cs[] <= 0. || cs[] >= 1.) {
64
65 /**
66 We need to make sure that the fine cells face fractions match
67 those of their neighbours. */
68
69 for (int j = 0; j <= 1; j++)
70 for (int k = 0; k <= 1; k++)
71 fine(fs.x,1,j,k) = cs[];
72 for (int i = 0; i <= 1; i++)
73 if (!is_refined(neighbor(2*i-1)) && neighbor(2*i-1).neighbors &&
74 (is_local(cell) || is_local(neighbor(2*i-1))))
75 for (int j = 0; j <= 1; j++)
76 for (int k = 0; k <= 1; k++)
77 fine(fs.x,2*i,j,k) = fs.x[i];
78 }
79 else {
80
81 /**
82 If the cell contains the embedded boundary, we reconstruct the
83 boundary using VOF linear reconstruction and a normal estimated
84 from the surface fractions. */
85
87 double alpha = plane_alpha (cs[], n);
88
89 /**
90 We need to reconstruct the face fractions *fs* for the fine cells.
91
92 For the fine face fractions contained within the coarse cell,
93 we compute the intersections directly using the VOF
94 reconstruction. */
95
96#if dimension == 2
97
98 /**
99 In 2D, we obtain the face fractions by taking into
100 account the orientation of the normal. */
101
102 if (2.*fabs(alpha) < fabs(n.y)) {
103 double yc = alpha/n.y;
104 int i = yc > 0.;
105 fine(fs.x,1,1 - i) = n.y < 0. ? 1. - i : i;
106 fine(fs.x,1,i) = n.y < 0. ? i - 2.*yc : 1. - i + 2.*yc;
107 }
108 else
109 fine(fs.x,1,0) = fine(fs.x,1,1) = alpha > 0.;
110
111#else // dimension == 3
112
113 /**
114 in 3D, we use the 2D projection of the reconstruction. */
115
116 for (int j = 0; j <= 1; j++)
117 for (int k = 0; k <= 1; k++)
118 if (!fine(cs,0,j,k) || !fine(cs,1,j,k))
119 fine(fs.x,1,j,k) = 0.;
120 else {
121 static const coord a = {0.,0.,0.}, b = {.5,.5,.5};
122 coord nc;
123 nc.x = 0., nc.y = (2.*j - 1.)*n.y, nc.z = (2.*k - 1.)*n.z;
124 fine(fs.x,1,j,k) = rectangle_fraction (nc, alpha, a, b);
125 }
126
127#endif // dimension == 3
128
129 /**
130 For the fine face fractions coincident with the faces of the
131 coarse cell, we obtain the intersection position from the
132 coarse cell face fraction. */
133
134 for (int i = 0; i <= 1; i++)
135 if (neighbor(2*i-1).neighbors &&
136 (is_local(cell) || is_local(neighbor(2*i-1)))) {
137 if (!is_refined(neighbor(2*i-1))) {
138 if (fs.x[i] <= 0. || fs.x[i] >= 1.)
139 for (int j = 0; j <= 1; j++)
140 for (int k = 0; k <= 1; k++)
141 fine(fs.x,2*i,j,k) = fs.x[i];
142 else {
143#if dimension == 2
144
145 /**
146 In 2D the orientation is obtained by looking at the values
147 of face fractions in the transverse direction. */
148
149 double a = fs.y[0,1] <= 0. || fs.y[2*i-1,1] <= 0. ||
150 fs.y[] >= 1. || fs.y[2*i-1] >= 1.;
151 if ((2.*a - 1)*(fs.x[i] - 0.5) > 0.) {
152 fine(fs.x,2*i,0) = a;
153 fine(fs.x,2*i,1) = 2.*fs.x[i] - a;
154 }
155 else {
156 fine(fs.x,2*i,0) = 2.*fs.x[i] + a - 1.;
157 fine(fs.x,2*i,1) = 1. - a;
158 }
159
160#else // dimension == 3
161
162 /**
163 In 3D we reconstruct the face fraction from the projection
164 of the cell interface reconstruction, as above. */
165
166 for (int j = 0; j <= 1; j++)
167 for (int k = 0; k <= 1; k++) {
168 static const coord a = {0.,0.,0.}, b = {.5,.5,.5};
169 coord nc;
170 nc.x = 0., nc.y = (2.*j - 1.)*n.y, nc.z = (2.*k - 1.)*n.z;
171 fine(fs.x,2*i,j,k) =
172 rectangle_fraction (nc, alpha - n.x*(2.*i - 1.)/2., a, b);
173 }
174
175#endif // dimension == 3
176 }
177 }
178
179 /**
180 The face fractions of empty children cells must be zero. */
181
182 for (int j = 0; j <= 1; j++)
183 #if dimension > 2
184 for (int k = 0; k <= 1; k++)
185 #endif
186 if (fine(fs.x,2*i,j,k) && !fine(cs,i,j,k))
187 fine(fs.x,2*i,j,k) = 0.;
188 }
189 }
190}
191
192/**
193## Restriction of cell-centered fields
194
195We now define restriction and prolongation functions for cell-centered
196fields. The goal is to define second-order operators which do not use
197any values from cells entirely contained within the embedded boundary
198(for which *cs = 0*).
199
200When restricting it is unfortunately not always possible to obtain a
201second-order interpolation. This happens when the parent cell does not
202contain enough child cells not entirely contained within the embedded
203boundary. In these cases, some external information (i.e. a boundary
204gradient condition) is required to be able to maintain second-order
205accuracy. This information can be passed by defining the
206*embed_gradient()* function of the field being restricted. */
207
210}
211
213{
214 // 0 children
215 if (!cs[]) {
216 s[] = 0.;
217 return;
218 }
219
220 /**
221 We first try to interpolate "diagonally". If enough child cells are
222 defined (i.e. have non-zero embedded fractions), we return the
223 corresponding value. */
224
225 double val = 0., nv = 0.;
226 for (int i = 0; i <= 1; i++)
227#if dimension > 2
228 for (int j = 0; j <= 1; j++)
229#endif
230 if (fine(cs,0,i,j) && fine(cs,1,!i,!j))
231 val += (fine(s,0,i,j) + fine(s,1,!i,!j))/2., nv++;
232 if (nv > 0.) {
233 s[] = val/nv;
234 return;
235 }
236
237 /**
238 Otherwise, we use the average of the child cells which are defined
239 (there is at least one). */
240
241 coord p = {0.,0.,0.};
242 for (int _c = 0; _c < 4; _c++) /* foreach_child */
243 if (cs[])
244 p.x += x, p.y += y, p.z += z, val += s[], nv++;
245 assert (nv > 0.);
246 s[] = val/nv;
247
248 /**
249 If the gradient is defined and if the variable is not using
250 homogeneous boundary conditions, we improve the interpolation using
251 this information. */
252
253 if (s.embed_gradient && s.boundary[0] != s.boundary_homogeneous[0]) {
254 coord o = {x,y,z}, g;
255 s.embed_gradient (point, s, &g);
256 for (int _d = 0; _d < dimension; _d++)
257 s[] += (o.x - p.x/nv)*g.x;
258 }
259}
260
261/**
262## Refinement/prolongation of cell-centered fields
263
264For refinement, we use either bilinear interpolation, if the required
265four coarse cell values are defined or trilinear interpolation if only
266three coarse cell values are defined. If less that three coarse cell
267values are defined ("pathological cases" below), we try to estimate
268gradients in each direction and add the corresponding correction. */
269
271{
272 for (int _c = 0; _c < 4; _c++) /* foreach_child */ {
273 if (!cs[])
274 s[] = 0.;
275 else {
276 assert (coarse(cs));
277 int i = (child.x + 1)/2, j = (child.y + 1)/2;
278#if dimension == 2
279 if (coarse(fs.x,i) && coarse(fs.y,0,j) &&
280 (coarse(cs) == 1. || coarse(cs,child.x) == 1. ||
281 coarse(cs,0,child.y) == 1. || coarse(cs,child.x,child.y) == 1.)) {
282 assert (coarse(cs,child.x) && coarse(cs,0,child.y));
283 if (coarse(fs.x,i,child.y) && coarse(fs.y,child.x,j)) {
284 // bilinear interpolation
285 assert (coarse(cs,child.x,child.y));
286 s[] = (9.*coarse(s) +
287 3.*(coarse(s,child.x) + coarse(s,0,child.y)) +
288 coarse(s,child.x,child.y))/16.;
289 }
290 else
291 // triangular interpolation
292 s[] = (2.*coarse(s) + coarse(s,child.x) + coarse(s,0,child.y))/4.;
293 }
294 else if (coarse(cs,child.x,child.y) &&
295 ((coarse(fs.x,i) && coarse(fs.y,child.x,j)) ||
296 (coarse(fs.y,0,j) && coarse(fs.x,i,child.y)))) {
297 // diagonal interpolation
298 s[] = (3.*coarse(s) + coarse(s,child.x,child.y))/4.;
299 }
300#else // dimension == 3
301 int k = (child.z + 1)/2;
302 if (coarse(fs.x,i) > 0.25 && coarse(fs.y,0,j) > 0.25 &&
303 coarse(fs.z,0,0,k) > 0.25 &&
304 (coarse(cs) == 1. || coarse(cs,child.x) == 1. ||
305 coarse(cs,0,child.y) == 1. || coarse(cs,child.x,child.y) == 1. ||
306 coarse(cs,0,0,child.z) == 1. || coarse(cs,child.x,0,child.z) == 1. ||
307 coarse(cs,0,child.y,child.z) == 1. ||
308 coarse(cs,child.x,child.y,child.z) == 1.)) {
309 assert (coarse(cs,child.x) && coarse(cs,0,child.y) &&
310 coarse(cs,0,0,child.z));
311 if (coarse(fs.x,i,child.y) && coarse(fs.y,child.x,j) &&
312 coarse(fs.z,child.x,child.y,k) &&
313 coarse(fs.z,child.x,0,k) && coarse(fs.z,0,child.y,k)) {
314 assert (coarse(cs,child.x,child.y) && coarse(cs,child.x,0,child.z) &&
315 coarse(cs,0,child.y,child.z) &&
316 coarse(cs,child.x,child.y,child.z));
317 // bilinear interpolation
318 s[] = (27.*coarse(s) +
319 9.*(coarse(s,child.x) + coarse(s,0,child.y) +
320 coarse(s,0,0,child.z)) +
321 3.*(coarse(s,child.x,child.y) + coarse(s,child.x,0,child.z) +
322 coarse(s,0,child.y,child.z)) +
323 coarse(s,child.x,child.y,child.z))/64.;
324 }
325 else
326 // tetrahedral interpolation
327 s[] = (coarse(s) + coarse(s,child.x) + coarse(s,0,child.y) +
328 coarse(s,0,0,child.z))/4.;
329 }
330 else if (coarse(cs,child.x,child.y,child.z) &&
331 ((coarse(fs.z,child.x,child.y,k) &&
332 ((coarse(fs.x,i) && coarse(fs.y,child.x,j)) ||
333 (coarse(fs.y,0,j) && coarse(fs.x,i,child.y))))
334 ||
335 (coarse(fs.z,0,0,k) &&
336 ((coarse(fs.x,i,0,child.z) && coarse(fs.y,child.x,j,child.z)) ||
337 (coarse(fs.y,0,j,child.z) && coarse(fs.x,i,child.y,child.z))))
338 ||
339 (coarse(fs.z,child.x,0,k) &&
340 coarse(fs.x,i) && coarse(fs.y,child.x,j,child.z))
341 ||
342 (coarse(fs.z,0,child.y,k) &&
343 coarse(fs.y,0,j) && coarse(fs.x,i,child.y,child.z))
344 ))
345 // diagonal interpolation
346 s[] = (3.*coarse(s) + coarse(s,child.x,child.y,child.z))/4.;
347#endif // dimension == 3
348 else {
349 // Pathological cases, use 1D gradients.
350 s[] = coarse(s);
351 for (int _d = 0; _d < dimension; _d++) {
352 if (coarse(fs.x,(child.x + 1)/2) && coarse(cs,child.x))
353 s[] += (coarse(s,child.x) - coarse(s))/4.;
354 else if (coarse(fs.x,(- child.x + 1)/2) && coarse(cs,- child.x))
355 s[] -= (coarse(s,- child.x) - coarse(s))/4.;
356 }
357 }
358 }
359 }
360}
361
362/**
363## Refinement/prolongation of face-centered velocity
364
365This function is modelled on
366[*refine_face_x()*](/src/grid/tree-common.h#refine_face_x) and is
367typically used to refine the values of the face-centered velocity
368field *uf*. It uses linear interpolation, taking into account the
369weighting by the embedded fractions *fs*. */
370
371for (int _d = 0; _d < dimension; _d++)
373{
374 vector v = s.v;
375 for (int i = 0; i <= 1; i++)
376 if (neighbor(2*i - 1).neighbors &&
377 (is_local(cell) || is_local(neighbor(2*i - 1)))) {
378 double g1 = fs.x[i] >= 1. && fs.x[i,+1] && fs.x[i,-1] ?
379 (v.x[i,+1]/fs.x[i,+1] - v.x[i,-1]/fs.x[i,-1])/8. : 0.;
380 double g2 = fs.x[i] >= 1. && fs.x[i,0,+1] && fs.x[i,0,-1] ?
381 (v.x[i,0,+1]/fs.x[i,0,+1] - v.x[i,0,-1]/fs.x[i,0,-1])/8. : 0.;
382 for (int j = 0; j <= 1; j++)
383 for (int k = 0; k <= 1; k++)
384 fine(v.x,2*i,j,k) = fs.x[i] ?
385 fine(fs.x,2*i,j,k)*(v.x[i]/fs.x[i] +
386 (2*j - 1)*g1 + (2*k - 1)*g2) : 0.;
387 }
389 double g1 = (fs.x[0,+1] + fs.x[1,+1]) && (fs.x[0,-1] + fs.x[1,-1]) ?
390 ((v.x[0,+1] + v.x[1,+1])/(fs.x[0,+1] + fs.x[1,+1]) -
391 (v.x[0,-1] + v.x[1,-1])/(fs.x[0,-1] + fs.x[1,-1]))/8. : 0.;
392 double g2 = (fs.x[1,0,+1] + fs.x[0,0,+1]) && (fs.x[1,0,-1] + fs.x[0,0,-1]) ?
393 ((v.x[0,0,+1] + v.x[1,0,+1])/(fs.x[1,0,+1] + fs.x[0,0,+1]) -
394 (v.x[0,0,-1] + v.x[1,0,-1])/(fs.x[1,0,-1] + fs.x[0,0,-1]))/8. : 0.;
395 for (int j = 0; j <= 1; j++)
396 for (int k = 0; k <= 1; k++)
397 fine(v.x,1,j,k) = fs.x[] + fs.x[1] ?
398 fine(fs.x,1,j,k)*((v.x[] + v.x[1])/(fs.x[] + fs.x[1]) +
399 (2*j - 1)*g1 + (2*k - 1)*g2) : 0.;
400 }
401}
vector g[]
We store the combined pressure gradient and acceleration field in g*.
Definition all-mach.h:65
const vector a
Definition all-mach.h:59
if TRASH define &&NewPid *& val(newpid, 0, 0, 0)) -> pid > 0) @else @ define is_newpid()(((NewPid *)&val(newpid, 0, 0, 0)) ->pid > 0) @endif Array *linear_tree(size_t size, scalar newpid)
Definition balance.h:13
#define dimension
Definition bitree.h:3
define k
define neighbor(o, p, q)((Point)
Definition cartesian.h:146
int y
Definition common.h:76
int x
Definition common.h:76
int z
Definition common.h:76
#define p
Definition tree.h:320
#define assert(a)
Definition config.h:107
Point point
Definition conserving.h:86
coord o
Definition curvature.h:672
else return n
Definition curvature.h:101
double alpha
Definition embed-tree.h:87
static void embed_fraction_refine(Point point, scalar cs)
Definition embed-tree.h:16
else fine(fs.x, 1, 0)
attribute
Definition embed-tree.h:208
static void restriction_embed_linear(Point point, scalar s)
Definition embed-tree.h:212
static void refine_embed_linear(Point point, scalar s)
Definition embed-tree.h:270
scalar s
Definition embed-tree.h:56
scalar cs[]
The volume and area fractions are stored in these fields.
Definition embed.h:21
static coord embed_gradient(Point point, vector u, coord p, coord n)
Definition embed.h:473
*cs[i, 0, 0] a *[i -1, 0, 0] j
Definition embed.h:88
double v[2]
Definition embed.h:383
vector fs[]
Definition embed.h:22
scalar int i
Definition embed.h:74
coord facet_normal(Point point, scalar c, vector s)
Definition fractions.h:426
double rectangle_fraction(coord n, double alpha, coord a, coord b)
VOF algorithms require the computation of volume fractions on (rectangular) parts of the initial squa...
Definition geometry.h:269
#define plane_alpha
Definition geometry.h:133
size *double * b
define is_active() cell(true) @define is_leaf(cell)(point.level
Definition linear.h:21
Definition common.h:78
double x
Definition common.h:79
scalar x
Definition common.h:47
scalar y
Definition common.h:49
define n n define coarse(a, k, p, n)((double *)(PARENT(k
#define is_refined(cell)
Definition tree.h:410
long nc
Definition utils.h:17
double cc
Definition vof.h:60