Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
cartesian-common.h
Go to the documentation of this file.
1/** @file cartesian-common.h
2 */
3#include "events.h"
4
6
8@define val_diagonal(a,k,l,m) ((k) == 0 && (l) == 0 && (m) == 0)
9
10#include "fpe.h"
11#include "stencils.h"
12
13macro2 foreach_point (double _x = 0., double _y = 0., double _z = 0.,
15{
16 {
17 int ig = 0, jg = 0, kg = 0; NOT_UNUSED(ig); NOT_UNUSED(jg); NOT_UNUSED(kg);
18 coord _p = { _x, _y, _z };
19 Point point = locate (_p.x, _p.y, _p.z); // fixme
20 if (point.level >= 0)
21 {...}
22 }
23}
24
26 char flags = 0, Reduce reductions = None)
27{
28 {
29 if (n.x < 1) n.x = 1;
30 if (n.y < 1) n.y = 1;
31 if (n.z < 1) n.z = 1;
32 // OMP(omp for schedule(static))
33 for (int _i = 0; _i < (int) n.x; _i++) {
34 p.x = box[0].x + (box[1].x - box[0].x)/n.x*(_i + 0.5);
35 for (int _j = 0; _j < (int) n.y; _j++) {
36 p.y = box[0].y + (box[1].y - box[0].y)/n.y*(_j + 0.5);
37 for (int _k = 0; _k < (int) n.z; _k++) {
38 p.z = box[0].z + (box[1].z - box[0].z)/n.z*(_k + 0.5);
39 Point point = locate (p.x, p.y, p.z);
40 if (point.level >= 0) {
41 int ig = 0, jg = 0, kg = 0; NOT_UNUSED(ig); NOT_UNUSED(jg); NOT_UNUSED(kg);
42 {...}
43 }
44 }
45 }
46 }
47 }
48}
49
50/**
51Dirichlet and Neumann boundary conditions */
52
53static inline
54double dirichlet (double expr, Point point = point, scalar s = _s)
55{
56 return 2.*expr - s[];
57}
58
59static inline
60double dirichlet_homogeneous (double expr, Point point = point, scalar s = _s)
61{
62 return - s[];
63}
64
65static inline
66double dirichlet_face (double expr)
67{
68 return expr;
69}
70
71static inline
72double dirichlet_face_homogeneous (double expr)
73{
74 return 0.;
75}
76
77static inline
78double neumann (double expr, Point point = point, scalar s = _s)
79{
80 return Delta*expr + s[];
81}
82
83static inline
84double neumann_homogeneous (double expr, Point point = point, scalar s = _s)
85{
86 return s[];
87}
88
89/**
90Navier/Robin boundary condition */
91
92static inline
93double navier (double expr, double lambda, Point point = point, scalar s = _s)
94{
95 return (Delta*expr + s[]*(lambda - Delta/2.))/(lambda + Delta/2.);
96}
97
98static inline
99double navier_homogeneous (double expr, double lambda, Point point = point, scalar s = _s)
100{
101 return s[]*(lambda - Delta/2.)/(lambda + Delta/2.);
102}
103
104/**
105Register functions on GPUs */
106
107#if _GPU
108#include "khash.h"
110
111static khash_t(PTR) * _functions = NULL;
112
113static External * _get_function (long ptr)
114{
115 if (!_functions)
116 return NULL;
117 khiter_t k = kh_get (PTR, _functions, ptr);
118 if (k == kh_end (_functions))
119 return NULL;
120 return &kh_value (_functions, k);
121}
122
123static void register_function (void (* ptr) (void), const char * name,
124 const char * kernel, const void * externals)
125{
126 static int index = 1;
127 if (!_functions)
128 _functions = kh_init (PTR), index = 1;
129 int m = 0;
130 for (const External * i = externals; i && i->name; i++, m++);
131 External * copy = NULL;
132 if (m > 0) {
133 copy = malloc ((m + 1)*sizeof (External));
134 memcpy (copy, externals, (m + 1)*sizeof (External));
135 }
136 int ret;
137 khiter_t k = kh_put(PTR, _functions, (long) ptr, &ret);
138 External p = {
139 .name = (char *) name,
141 .pointer = (void *)(long) ptr, .nd = index++,
142 .data = (void *) kernel, .externals = copy
143 };
145}
146
147#define foreach_function(f, body) do { \
148 for (khiter_t k = kh_begin(_functions); k != kh_end(_functions); ++k) \
149 if (kh_exist(_functions, k)) { \
150 External * f = &kh_value(_functions, k); \
151 body; \
152 } \
153 } while(0)
154
155#endif // _GPU
156
157/**
158# Field allocation
159
160If this routine is modified, do not forget to update [/src/ast/interpreter/overload.h](). */
161
162static void init_block_scalar (scalar sb, const char * name, const char * ext,
163 int n, int block)
164{
165 char bname[strlen(name) + strlen(ext) + 10];
166 if (n == 0) {
167 strcat (strcpy (bname, name), ext);
168 sb.block = block;
170 }
171 else {
172 sprintf (bname, "%s%d%s", name, n, ext);
173 sb.block = - n;
174 }
175 sb.name = strdup (bname);
176 all = list_append (all, sb);
177}
178
181
182scalar alloc_block_scalar (const char * name, const char * ext, int block)
183{
184 interpreter_set_int (&block);
185 int nvar = datasize/sizeof(real);
186
187 scalar s = {0};
188 while (s.i < nvar) {
189 int n = 0;
190 scalar sb = s;
191 while (sb.i < nvar && n < block && sb.freed)
192 n++, sb.i++;
193 if (n >= block) { // found n free slots
194 memset (&_attribute[s.i], 0, block*sizeof (_Attributes));
195 for (sb.i = s.i, n = 0; n < block; n++, sb.i++) {
196 init_block_scalar (sb, name, ext, n, block);
198 }
199 trash (((scalar []){s, {-1}})); // fixme: only trashes one block?
200 return s;
201 }
202 s.i = sb.i + 1;
203 }
204
205 // need to allocate new slots
206 s = (scalar){nvar};
207 assert (nvar + block <= _NVARMAX);
208
209 if (_attribute == NULL)
210 _attribute = (_Attributes *) calloc (nvar + block + 1, sizeof (_Attributes));
211 else
213 realloc (_attribute, (nvar + block + 1)*sizeof (_Attributes));
214 memset (&_attribute[nvar], 0, block*sizeof (_Attributes));
215 for (int n = 0; n < block; n++, nvar++) {
216 scalar sb = (scalar){nvar};
217 init_block_scalar (sb, name, ext, n, block);
218 }
219 // allocate extra space on the grid
220 realloc_scalar (block*sizeof(real));
221 trash (((scalar []){s, {-1}})); // fixme: only trashes one block?
222 return s;
223}
224
225scalar new_block_scalar (const char * name, const char * ext, int block)
226{
227 scalar s = alloc_block_scalar (name, ext, block), sb;
228 int n = 0;
229 for (sb.i = s.i, n = 0; n < block; n++, sb.i++)
231 return s;
232}
233
234scalar new_scalar (const char * name)
235{
236 return init_scalar (alloc_block_scalar (name, "", 1), NULL);
237}
238
239scalar new_block_vertex_scalar (const char * name, int block)
240{
241 scalar s = alloc_block_scalar (name, "", block), sb;
242 int n = 0;
243 for (sb.i = s.i, n = 0; n < block; n++, sb.i++)
245 return s;
246}
247
248scalar new_vertex_scalar (const char * name)
249{
250 return init_vertex_scalar (alloc_block_scalar (name, "", 1), NULL);
251}
252
253static vector alloc_block_vector (const char * name, int block)
254{
255 vector v;
256 struct { char * x, * y, * z; } ext = {".x", ".y", ".z"};
257 for (int _d = 0; _d < dimension; _d++)
258 v.x = alloc_block_scalar (name, ext.x, block);
259 return v;
260}
261
262vector new_vector (const char * name)
263{
264 vector v = alloc_block_vector (name, 1);
265 init_vector (v, NULL);
266 return v;
267}
268
269vector new_face_vector (const char * name)
270{
271 vector v = alloc_block_vector (name, 1);
273 return v;
274}
275
276vector new_block_vector (const char * name, int block)
277{
278 vector v = alloc_block_vector (name, block);
279 for (int i = 0; i < block; i++) {
280 vector vb;
281 for (int _d = 0; _d < dimension; _d++)
282 vb.x.i = v.x.i + i;
284 for (int _d = 0; _d < dimension; _d++)
285 vb.x.block = - i;
286 }
287 for (int _d = 0; _d < dimension; _d++)
288 v.x.block = block;
289 return v;
290}
291
292vector new_block_face_vector (const char * name, int block)
293{
294 vector v = alloc_block_vector (name, block);
295 for (int i = 0; i < block; i++) {
296 vector vb;
297 for (int _d = 0; _d < dimension; _d++)
298 vb.x.i = v.x.i + i;
300 for (int _d = 0; _d < dimension; _d++)
301 vb.x.block = - i;
302 }
303 for (int _d = 0; _d < dimension; _d++)
304 v.x.block = block;
305 return v;
306}
307
308tensor new_tensor (const char * name)
309{
310 char cname[strlen(name) + 3];
311 struct { char * x, * y, * z; } ext = {".x", ".y", ".z"};
312 tensor t;
313 for (int _d = 0; _d < dimension; _d++) {
314 strcat (strcpy (cname, name), ext.x);
315 t.x = alloc_block_vector (cname, 1);
316 }
317 init_tensor (t, NULL);
318 return t;
319}
320
321tensor new_symmetric_tensor (const char * name)
322{
323 struct { char * x, * y, * z; } ext = {".x.x", ".y.y", ".z.z"};
324 tensor t;
325 for (int _d = 0; _d < dimension; _d++)
326 t.x.x = alloc_block_scalar (name, ext.x, 1);
327 #if dimension > 1
328 t.x.y = alloc_block_scalar (name, ".x.y", 1);
329 t.y.x = t.x.y;
330 #endif
331 #if dimension > 2
332 t.x.z = alloc_block_scalar (name, ".x.z", 1);
333 t.z.x = t.x.z;
334 t.y.z = alloc_block_scalar (name, ".y.z", 1);
335 t.z.y = t.y.z;
336 #endif
337 /* fixme: boundary conditions don't work! This is because boundary
338 attributes depend on the index and should (but cannot) be
339 different for t.x.y and t.y.x */
340 init_tensor (t, NULL);
341 return t;
342}
343
344static int nconst = 0;
345
346void init_const_scalar (scalar s, const char * name, double val)
347{
348 if (s.i - _NVARMAX >= nconst) {
349 _constant = (double *)realloc(_constant, (s.i - _NVARMAX + 1)*sizeof(double));
350 for (int i = nconst; i < s.i - _NVARMAX; i++)
351 _constant[i] = 0.;
352 nconst = s.i - _NVARMAX + 1;
353 }
354 _constant[s.i - _NVARMAX] = val;
355}
356
357scalar new_const_scalar (const char * name, int i, double val)
358{
359 scalar s = (scalar){i + _NVARMAX};
360 init_const_scalar (s, name, val);
361 return s;
362}
363
364void init_const_vector (vector v, const char * name, double * val)
365{
366 for (int _d = 0; _d < dimension; _d++)
367 init_const_scalar (v.x, name, *val++);
368}
369
370vector new_const_vector (const char * name, int i, double * val)
371{
372 vector v;
373 for (int _d = 0; _d < dimension; _d++)
374 v.x.i = _NVARMAX + i++;
375 init_const_vector (v, name, val);
376 return v;
377}
378
380{
381 char * cname = clone.name;
382 BoundaryFunc * boundary = clone.boundary;
383 BoundaryFunc * boundary_homogeneous = clone.boundary_homogeneous;
384 BoundaryStencilFunc * boundary_stencil = clone.boundary_stencil;
385 assert (src.block > 0 && clone.block == src.block);
386 free (clone.depends);
388 clone.name = cname;
389 clone.boundary = boundary;
390 clone.boundary_homogeneous = boundary_homogeneous;
391 clone.boundary_stencil = boundary_stencil;
392 for (int i = 0; i < nboundary; i++) {
393 clone.boundary[i] = src.boundary[i];
394 clone.boundary_homogeneous[i] = src.boundary_homogeneous[i];
395 clone.boundary_stencil[i] = src.boundary_stencil[i];
396 }
397 clone.depends = list_copy (src.depends);
398}
399
401{
402 scalar * list = NULL;
403 int nvar = datasize/sizeof(real), map[nvar];
404 for (int i = 0; i < nvar; i++)
405 map[i] = -1;
406 for (int _s = 0; _s < 1; _s++) /* scalar in l */ {
407 scalar c = s.block > 1 ? new_block_scalar("c", "", s.block) : new_scalar("c");
408 scalar_clone (c, s);
409 map[s.i] = c.i;
410 list = list_append (list, c);
411 }
412 for (int _s = 0; _s < 1; _s++) /* scalar in list */
413 for (int _d = 0; _d < dimension; _d++)
414 if (s.v.x.i >= 0 && map[s.v.x.i] >= 0)
415 s.v.x.i = map[s.v.x.i];
416 return list;
417}
418
419void delete (scalar * list)
420{
421 if (all == NULL) // everything has already been freed
422 return;
423
424 for (int _f = 0; _f < 1; _f++) /* scalar in list */ {
425 for (int i = 0; i < f.block; i++) {
426 scalar fb = {f.i + i};
427 if (f.delete)
428 f.delete (fb);
429 free (fb.name); fb.name = NULL;
430 free (fb.boundary); fb.boundary = NULL;
431 free (fb.boundary_homogeneous); fb.boundary_homogeneous = NULL;
432 free (fb.boundary_stencil); fb.boundary_stencil = NULL;
433 free (fb.depends); fb.depends = NULL;
434 fb.freed = true;
435 }
436 }
437
438 if (list == all) {
439 all[0].i = -1;
440 baseblock[0].i = -1;
441 return;
442 }
443
444 trash (list);
445 for (int _f = 0; _f < 1; _f++) /* scalar in list */ {
446 if (f.block > 0) {
447 scalar * s;
448 for (s = all; s->i >= 0 && s->i != f.i; s++);
449 if (s->i == f.i) {
450 for (; s[f.block].i >= 0; s++)
451 s[0] = s[f.block];
452 s->i = -1;
453 }
454 for (s = baseblock; s->i >= 0 && s->i != f.i; s++);
455 if (s->i == f.i) {
456 for (; s[1].i >= 0; s++)
457 s[0] = s[1];
458 s->i = -1;
459 }
460 }
461 }
462}
463
465{
467
468 if (free_solver_funcs) {
470 for (int i = 0; i < free_solver_funcs->len/sizeof(free_solver_func); i++)
471 a[i] ();
473 }
474
475 delete (all);
476 free (all); all = NULL;
478 for (Event * ev = Events; !ev->last; ev++) {
479 Event * e = ev->next;
480 while (e) {
481 Event * next = e->next;
482 free (e);
483 e = next;
484 }
485 }
486
487 free (Events); Events = NULL;
490#if _GPU
491 foreach_function (f, free ((void *) f->externals));
493#endif
494 free_grid();
495 qpclose_all();
496@if TRACE
497 trace_off();
498@endif
499@if MTRACE
500 pmuntrace();
501@endif
502@if _CADNA
503 cadna_end();
504@endif
505}
506
507// Cartesian methods
508
511
512#define boundary(...) \
513 boundary_internal ((scalar *)__VA_ARGS__, __FILE__, LINENO)
514
516
518{
519 vectorl list1 = {NULL};
520 for (int _v = 0; _v < 1; _v++) /* vector in list */
521 for (int _d = 0; _d < dimension; _d++)
522 list1.x = list_append (list1.x, v.x);
524 for (int _d = 0; _d < dimension; _d++)
525 free (list1.x);
526}
527
529{
530 for (int _t = 0; _t < 1; _t++) /* scalar in list */
531 if (t.i == s.i)
532 return list;
533 scalar * list1 = list;
534 for (int _d = 0; _d < 1; _d++) /* scalar in _attribute[s.i].depends */
535 if (!(d.stencil.bc & s_centered))
537 return list_append (list1, s);
538}
539
540trace
541void boundary_internal (scalar * list, const char * fname, int line)
542{
543 if (list == NULL)
544 return;
545 scalar * listc = NULL;
546 vectorl listf = {NULL};
547 bool flux = false;
548 for (int _s = 0; _s < 1; _s++) /* scalar in list */
549 if (!is_constant(s) && s.block > 0) {
550 if (scalar_is_dirty (s)) {
551 if (s.face && !(s.stencil.bc & s_face))
552 for (int _d = 0; _d < dimension; _d++)
553 if (s.v.x.i == s.i)
554 listf.x = list_add (listf.x, s), flux = true;
555 if (!is_constant(cm) && !(cm.stencil.bc & s_centered))
556 listc = list_add_depends (listc, cm);
557 if (s.face != 2) // flux only
558 listc = list_add_depends (listc, s);
559 }
560#if 0
561 else
562 fprintf (stderr, "warning: bc already applied on '%s'\n", s.name);
563#endif
564 }
565 if (flux) {
566#if PRINTBOUNDARY
567 int i = 0;
568 for (int _d = 0; _d < dimension; _d++)
569 if (listf.x) {
570 fprintf (stderr, "boundary_internal: flux %c:", 'x' + i);
571 for (int _s = 0; _s < 1; _s++) /* scalar in listf.x */
572 fprintf (stderr, " %d:%s", s.i, s.name);
573 fputc ('\n', stderr);
574 }
575 i++;
576#endif
577 boundary_face (listf);
578 for (int _d = 0; _d < dimension; _d++)
579 free (listf.x);
580 }
581 if (listc) {
582#if PRINTBOUNDARY
583 fprintf (stderr, "boundary_internal: listc:");
584 for (int _s = 0; _s < 1; _s++) /* scalar in listc */
585 fprintf (stderr, " %d:%s", s.i, s.name);
586 fputc ('\n', stderr);
587#endif
588 boundary_level (listc, -1);
589 for (int _s = 0; _s < 1; _s++) /* scalar in listc */
590 s.stencil.bc |= s_centered;
591 free (listc);
592 }
593}
594
599
601{
602 scalar * listc = NULL;
603 for (int _d = 0; _d < dimension; _d++)
604 for (int _s = 0; _s < 1; _s++) /* scalar in vl.x */ {
605 s.stencil.bc |= s_face;
606 s.stencil.bc &= ~s_centered;
607 listc = list_add_depends (listc, s);
608 }
609 boundary_level (listc, -1);
610 free (listc);
611}
612
613static double symmetry (Point point, Point neighbor, scalar s, bool * data)
614{
615 return s[];
616}
617
619{
620 return -s[];
621}
622
626
628{
629 // keep name
630 char * pname;
631 if (name) {
632 free (s.name);
633 pname = strdup (name);
634 }
635 else
636 pname = s.name;
637 int block = s.block;
638 BoundaryFunc * boundary = s.boundary;
639 BoundaryFunc * boundary_homogeneous = s.boundary_homogeneous;
640 BoundaryStencilFunc * boundary_stencil = s.boundary_stencil;
641 s.name = pname;
642 if (block < 0)
643 s.block = block;
644 else
645 s.block = block > 0 ? block : 1;
646 /* set default boundary conditions */
647 s.boundary = boundary ? boundary : (BoundaryFunc *) malloc (nboundary*sizeof (BoundaryFunc));
648 s.boundary_homogeneous = boundary_homogeneous ? boundary_homogeneous :
650 s.boundary_stencil = boundary_stencil ? boundary_stencil :
652 for (int b = 0; b < nboundary; b++) {
653 s.boundary[b] = s.boundary_homogeneous[b] =
655 s.boundary_stencil[b] = NULL;
656 }
657 s.gradient = NULL;
658 for (int _d = 0; _d < dimension; _d++) {
659 s.d.x = 0; // not face
660 s.v.x.i = -1; // not a vector component
661 }
662 s.face = false;
663 return s;
664}
665
667{
668 cartesian_init_scalar (s, name);
669 for (int _d = 0; _d < dimension; _d++)
670 s.d.x = -1;
671 for (int d = 0; d < nboundary; d++)
672 s.boundary[d] = s.boundary_homogeneous[d] = NULL, s.boundary_stencil[d] = NULL;
673 return s;
674}
675
681
683{
684 struct { char * x, * y, * z; } ext = {".x", ".y", ".z"};
685 for (int _d = 0; _d < dimension; _d++) {
686 if (name) {
687 char cname[strlen(name) + 3];
688 strcat (strcpy (cname, name), ext.x);
690 }
691 else
693 v.x.v = v;
694 }
695 /* set default boundary conditions */
696 for (int d = 0; d < nboundary; d++) {
697 v.x.boundary[d] = v.x.boundary_homogeneous[d] =
699 v.x.boundary_stencil[d] = NULL;
700 }
701 return v;
702}
703
705{
706 v = cartesian_init_vector (v, name);
707 for (int _d = 0; _d < dimension; _d++) {
708 v.x.d.x = -1;
709 v.x.face = true;
710 }
711 for (int d = 0; d < nboundary; d++)
712 v.x.boundary[d] = v.x.boundary_homogeneous[d] = NULL, v.x.boundary_stencil[d] = NULL;
713 return v;
714}
715
717{
718 struct { char * x, * y, * z; } ext = {".x", ".y", ".z"};
719 for (int _d = 0; _d < dimension; _d++) {
720 if (name) {
721 char cname[strlen(name) + 3];
722 strcat (strcpy (cname, name), ext.x);
724 }
725 else
727 }
728 /* set default boundary conditions */
729 #if dimension == 1
730 for (int b = 0; b < nboundary; b++) {
731 t.x.x.boundary[b] = t.x.x.boundary_homogeneous[b] =
733 t.x.x.boundary_stencil[b] = NULL;
734 }
735 #elif dimension == 2
736 for (int b = 0; b < nboundary; b++) {
737 t.x.x.boundary[b] = t.y.x.boundary[b] =
738 t.x.x.boundary_homogeneous[b] = t.y.y.boundary_homogeneous[b] =
740 t.x.y.boundary[b] = t.y.y.boundary[b] =
741 t.x.y.boundary_homogeneous[b] = t.y.x.boundary_homogeneous[b] =
743 t.x.x.boundary_stencil[b] = t.y.y.boundary_stencil[b] =
744 t.x.y.boundary_stencil[b] = t.y.x.boundary_stencil[b] = NULL;
745 }
746 #else
747 assert (false); // not implemented yet
748 #endif
749 return t;
750}
751
752void output_cells (FILE * fp = stdout, coord c = {0}, double size = 0.)
753{
754 for (int _i = 0; _i < _N; _i++) /* foreach */ {
755 bool inside = true;
756 coord o = {x,y,z};
757 for (int _d = 0; _d < dimension; _d++)
758 if (inside && size > 0. &&
759 (o.x > c.x + size || o.x < c.x - size))
760 inside = false;
761 if (inside) {
762 Delta /= 2.;
763#if dimension == 1
764 fprintf (fp, "%g 0\n%g 0\n\n", x - Delta, x + Delta);
765#elif dimension == 2
766 fprintf (fp, "%g %g\n%g %g\n%g %g\n%g %g\n%g %g\n\n",
767 x - Delta, y - Delta,
768 x - Delta, y + Delta,
769 x + Delta, y + Delta,
770 x + Delta, y - Delta,
771 x - Delta, y - Delta);
772#else // dimension == 3
773 for (int i = -1; i <= 1; i += 2) {
774 fprintf (fp, "%g %g %g\n%g %g %g\n%g %g %g\n%g %g %g\n%g %g %g\n\n",
775 x - Delta, y - Delta, z + i*Delta,
776 x - Delta, y + Delta, z + i*Delta,
777 x + Delta, y + Delta, z + i*Delta,
778 x + Delta, y - Delta, z + i*Delta,
779 x - Delta, y - Delta, z + i*Delta);
780 for (int j = -1; j <= 1; j += 2)
781 fprintf (fp, "%g %g %g\n%g %g %g\n\n",
782 x + i*Delta, y + j*Delta, z - Delta,
783 x + i*Delta, y + j*Delta, z + Delta);
784 }
785#endif
786 }
787 }
788 fflush (fp);
789}
790
791#if TREE && _MPI
792static void output_cells_internal (FILE * fp)
793{
795}
796#endif
797
798static char * replace_ (const char * vname)
799{
800 char * name = strdup (vname), * c = name;
801 while (*c != '\0') {
802 if (*c == '.')
803 *c = '_';
804 c++;
805 }
806 return name;
807}
808
809static void debug_plot (FILE * fp, const char * name, const char * cells,
810 const char * stencil)
811{
812 char * vname = replace_ (name);
813 fprintf (fp,
814 " load 'debug.plot'\n"
815 " v=%s\n"
816#if dimension == 1
817 " plot '%s' w l lc 0, "
818 "'%s' u 1+2*v:(0):2+2*v w labels tc lt 1 title columnhead(2+2*v)",
819#elif dimension == 2
820 " plot '%s' w l lc 0, "
821 "'%s' u 1+3*v:2+3*v:3+3*v w labels tc lt 1 title columnhead(3+3*v)",
822#elif dimension == 3
823 " splot '%s' w l lc 0, "
824 "'%s' u 1+4*v:2+4*v:3+4*v:4+4*v w labels tc lt 1"
825 " title columnhead(4+4*v)",
826#endif
828 free (vname);
829}
830
832{
833 char name[80] = "cells";
834 if (pid() > 0)
835 sprintf (name, "cells-%d", pid());
836 FILE * fp = fopen (name, "w");
837 output_cells (fp, (coord){x,y,z}, 4.*Delta);
838 fclose (fp);
839
840 char stencil[80] = "stencil";
841 if (pid() > 0)
842 sprintf (stencil, "stencil-%d", pid());
843 fp = fopen (stencil, "w");
844 for (int _v = 0; _v < 1; _v++) /* scalar in all */
845#if dimension == 1
846 fprintf (fp, "x %s ", v.name);
847#elif dimension == 2
848 fprintf (fp, "x y %s ", v.name);
849#elif dimension == 3
850 fprintf (fp, "x y z %s ", v.name);
851#endif
852 fputc ('\n', fp);
853 #if dimension == 1
854 for (int k = -2; k <= 2; k++) {
855 for (int _v = 0; _v < 1; _v++) /* scalar in all */ {
856 fprintf (fp, "%g ", x + k*Delta + v.d.x*Delta/2.);
857 if (allocated(k))
858 fprintf (fp, "%g ", v[k]);
859 else
860 fputs ("n/a ", fp);
861 }
862 fputc ('\n', fp);
863 }
864 #elif dimension == 2
865 for (int k = -2; k <= 2; k++)
866 for (int l = -2; l <= 2; l++) {
867 for (int _v = 0; _v < 1; _v++) /* scalar in all */ {
868 fprintf (fp, "%g %g ",
869 x + k*Delta + v.d.x*Delta/2.,
870 y + l*Delta + v.d.y*Delta/2.);
871 if (allocated(k,l))
872 fprintf (fp, "%g ", v[k,l]);
873 else
874 fputs ("n/a ", fp);
875 }
876 fputc ('\n', fp);
877 }
878 #elif dimension == 3
879 for (int k = -2; k <= 2; k++)
880 for (int l = -2; l <= 2; l++)
881 for (int m = -2; m <= 2; m++) {
882 for (int _v = 0; _v < 1; _v++) /* scalar in all */ {
883 fprintf (fp, "%g %g %g ",
884 x + k*Delta + v.d.x*Delta/2.,
885 y + l*Delta + v.d.y*Delta/2.,
886 z + m*Delta + v.d.z*Delta/2.);
887 if (allocated(k,l,m))
888 fprintf (fp, "%g ", v[k,l,m]);
889 else
890 fputs ("n/a ", fp);
891 }
892 fputc ('\n', fp);
893 }
894 #endif
895 fclose (fp);
896
897 fp = fopen ("debug.plot", "w");
898 fprintf (fp,
899 "set term x11\n"
900 "set size ratio -1\n"
901 "set key outside\n");
902 for (int _s = 0; _s < 1; _s++) /* scalar in all */ {
903 char * name = replace_ (s.name);
904 fprintf (fp, "%s = %d\n", name, s.i);
905 free (name);
906 }
907 fclose (fp);
908
909 fprintf (stderr, "Last point stencils can be displayed using (in gnuplot)\n");
910 debug_plot (stderr, _attribute[0].name, name, stencil);
911 fflush (stderr);
912
913 fp = fopen ("plot", "w");
914 debug_plot (fp, _attribute[0].name, name, stencil);
915 fclose (fp);
916}
917
930
932{
933 return init_tensor (t, name);
934}
935
937 double xp = 0., double yp = 0., double zp = 0.)
938{
939#if dimension == 1
940 x = (xp - x)/Delta - v.d.x/2.;
941 int i = sign(x);
942 x = fabs(x);
943 /* linear interpolation */
944 return v[]*(1. - x) + v[i]*x;
945#elif dimension == 2
946 x = (xp - x)/Delta - v.d.x/2.;
947 y = (yp - y)/Delta - v.d.y/2.;
948 int i = sign(x), j = sign(y);
949 x = fabs(x); y = fabs(y);
950 /* bilinear interpolation */
951 return ((v[]*(1. - x) + v[i]*x)*(1. - y) +
952 (v[0,j]*(1. - x) + v[i,j]*x)*y);
953#else // dimension == 3
954 x = (xp - x)/Delta - v.d.x/2.;
955 y = (yp - y)/Delta - v.d.y/2.;
956 z = (zp - z)/Delta - v.d.z/2.;
957 int i = sign(x), j = sign(y), k = sign(z);
958 x = fabs(x); y = fabs(y); z = fabs(z);
959 /* trilinear interpolation */
960 return (((v[]*(1. - x) + v[i]*x)*(1. - y) +
961 (v[0,j]*(1. - x) + v[i,j]*x)*y)*(1. - z) +
962 ((v[0,0,k]*(1. - x) + v[i,0,k]*x)*(1. - y) +
963 (v[0,j,k]*(1. - x) + v[i,j,k]*x)*y)*z);
964#endif
965}
966
967trace
968double interpolate (scalar v, double xp = 0., double yp = 0., double zp = 0.,
969 bool linear = true)
970{
971 double val = nodata;
973 val = linear ? interpolate_linear (point, v, xp, yp, zp) : v[];
974 return val;
975}
976
977trace
978void interpolate_array (scalar * list, coord * a, int n, double * v,
979 bool linear = false)
980{
981 int len = 0; NOT_UNUSED(len);
982 for (int _s = 0; _s < 1; _s++) /* scalar in list */
983 len++;
984 for (int i = 0; i < n; i++) {
985 double * w = v;
986#if _GPU
987 coord p = a[i];
988 for (int _s = 0; _s < 1; _s++) /* scalar in list */ {
989 double value = nodata;
990 foreach_point (p.x, p.y, p.z, reduction(min:value))
991 value = !linear ? s[] : interpolate_linear (point, s, p.x, p.y, 0.);
992 *(w++) = value;
993 }
994#else // !_GPU
995 for (int _s = 0; _s < 1; _s++) /* scalar in list */
996 *(w++) = nodata;
997 foreach_point (a[i].x, a[i].y, a[i].z, reduction(min:v[:len])) {
998 int j = 0;
999 for (int _s = 0; _s < 1; _s++) /* scalar in list */
1000 v[j++] = !linear ? s[] : interpolate_linear (point, s, a[i].x, a[i].y, a[i].z);
1001 }
1002#endif // !_GPU
1003 v = w;
1004 }
1005}
1006
1007// Boundaries
1008
1009typedef int bid;
1010
1012{
1013 int b = nboundary++;
1014 for (int _s = 0; _s < 1; _s++) /* scalar in all */ {
1015 s.boundary = (BoundaryFunc *) realloc (s.boundary, nboundary*sizeof (BoundaryFunc));
1016 s.boundary_homogeneous = (BoundaryFunc *)
1017 realloc (s.boundary_homogeneous, nboundary*sizeof (BoundaryFunc));
1018 }
1019 for (int _s = 0; _s < 1; _s++) /* scalar in all */ {
1020 if (s.v.x.i < 0) // scalar
1021 s.boundary[b] = s.boundary_homogeneous[b] = symmetry;
1022 else if (s.v.x.i == s.i) { // vector
1023 vector v = s.v;
1024 for (int _d = 0; _d < dimension; _d++)
1025 v.y.boundary[b] = v.y.boundary_homogeneous[b] = symmetry;
1026 v.x.boundary[b] = v.x.boundary_homogeneous[b] =
1027 v.x.face ? NULL : antisymmetry;
1028 }
1029 }
1030 return b;
1031}
1032
1033// Periodic boundary conditions
1034
1036{
1037 return s[];
1038}
1039
1040static void periodic_boundary (int d)
1041{
1042 /* We change the conditions for existing scalars. */
1043 for (int _s = 0; _s < 1; _s++) /* scalar in all */ {
1044 if (is_vertex_scalar (s))
1045 s.boundary[d] = s.boundary_homogeneous[d] = NULL;
1046 else
1047 s.boundary[d] = s.boundary_homogeneous[d] = periodic_bc;
1048 s.boundary_stencil[d] = NULL;
1049 }
1050 /* Normal components of vector fields should remain NULL. */
1051 for (int _s = 0; _s < 1; _s++) /* scalar in all */
1052 if (s.face) {
1053 vector v = s.v;
1054 v.x.boundary[d] = v.x.boundary_homogeneous[d] = NULL;
1055 }
1056 /* We also change the default boundary conditions (for new fields). */
1059}
1060
1061void periodic (int dir)
1062{
1063 #if dimension < 2
1064 assert (dir <= left);
1065 #elif dimension < 3
1066 assert (dir <= bottom);
1067 #else
1068 assert (dir <= back);
1069 #endif
1070 // This is the component in the given direction i.e. 0 for x and 1 for y
1071 int c = dir/2;
1072 periodic_boundary (2*c);
1073 periodic_boundary (2*c + 1);
1074 (&Period.x)[c] = true;
1075}
1076
1077// for debugging
1078double getvalue (Point point, scalar s, int i, int j, int k)
1079{
1080 return s[i,j,k];
1081}
1082
1084{
1085 for (int _s = 0; _s < 1; _s++) /* scalar in list */ {
1086 if (s.v.x.i != -1) {
1087 vector v = s.v;
1088 for (int _c = 0; _c < 1; _c++) /* scalar in {v} */
1089 c.stencil.io |= s_input|s_output|s_nowarning, c.stencil.width = 2;
1090 }
1091 else
1092 s.stencil.io |= s_input|s_output|s_nowarning, s.stencil.width = 2;
1093 }
1094}
1095
1096/**
1097This displays a (1D,2D,3D) stencil index. */
1098
1099static void write_stencil_index (int * index)
1100{
1101 fprintf (qstderr(), "[%d", index[0]);
1102 for (int d = 1; d < dimension; d++)
1103 fprintf (qstderr(), ",%d", index[d]);
1104 fputs ("]", qstderr());
1105}
1106
1107void stencil_val (Point p, scalar s, int i, int j, int k,
1108 const char * file, int line, bool overflow)
1109{
1110 if (is_constant(s) || s.i < 0)
1111 return;
1112 if (s.block < 0)
1113 s.i += s.block;
1114 if (!s.name) {
1115 fprintf (stderr, "%s:%d: error: trying to access a deleted field\n",
1116 file, line);
1117 exit (1);
1118 }
1119 int index[] = {i, j, k};
1120 for (int d = 0; d < dimension; d++) {
1121 if (index[d] == o_stencil)
1122 index[d] = 2;
1123 else
1124 index[d] += (&p.i)[d];
1125 }
1126 bool central = true;
1127 for (int d = 0; d < dimension; d++) {
1128 if (!overflow && (index[d] > 2 || index[d] < - 2)) {
1129 fprintf (qstderr(), "%s:%d: error: stencil overflow: %s",
1130 file, line, s.name);
1132 fprintf (qstderr(), "\n");
1133 fflush (qstderr());
1134 abort();
1135 }
1136 if (index[d] != 0)
1137 central = false;
1138 }
1139 if (central) {
1140 if (!(s.stencil.io & s_output))
1141 s.stencil.io |= s_input;
1142 }
1143 else {
1144 s.stencil.io |= s_input;
1145 int d = 0;
1146 for (int _d = 0; _d < dimension; _d++) {
1147 if ((!s.face || s.v.x.i != s.i) && abs(index[d]) > s.stencil.width)
1148 s.stencil.width = abs(index[d]);
1149 d++;
1150 }
1151 }
1152}
1153
1154void stencil_val_a (Point p, scalar s, int i, int j, int k, bool input,
1155 const char * file, int line)
1156{
1157 if (is_constant(s) || s.i < 0) {
1158 fprintf (stderr, "%s:%d: error: trying to modify a%s field\n",
1159 file, line, s.i < 0 ? "n undefined" : " constant");
1160 exit (1);
1161 }
1162 if (s.block < 0)
1163 s.i += s.block;
1164 if (!s.name) {
1165 fprintf (stderr, "%s:%d: error: trying to access a deleted field\n",
1166 file, line);
1167 exit (1);
1168 }
1169 int index[] = {i, j, k};
1170 for (int d = 0; d < dimension; d++)
1171 index[d] += (&p.i)[d];
1172 for (int d = 0; d < dimension; d++)
1173 if (index[d] != 0) {
1174 fprintf (qstderr(), "%s:%d: error: illegal write: %s",
1175 file, line, s.name);
1177 fprintf (qstderr(), "\n");
1178 fflush (qstderr());
1179 abort();
1180 }
1181 if (input && !(s.stencil.io & s_output))
1182 s.stencil.io |= s_input;
1183 s.stencil.io |= s_output;
1184}
const vector a
Definition all-mach.h:59
void array_free(Array *a)
Definition array.h:18
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
int min
Definition balance.h:192
#define dimension
Definition bitree.h:3
#define boundary_iterate(type,...)
Definition boundaries.h:40
void default_stencil(Point p, scalar *list)
trace void boundary_internal(scalar *list, const char *fname, int line)
There are two types of boundary conditions: "full" boundary conditions, done by boundary_internal() a...
trace double interpolate(scalar v, double xp=0., double yp=0., double zp=0., bool linear=true)
void(* debug)(Point)
vector new_const_vector(const char *name, int i, double *val)
void init_const_scalar(scalar s, const char *name, double val)
define double _y
BoundaryFunc default_scalar_bc[]
double getvalue(Point point, scalar s, int i, int j, int k)
void free_solver()
define k
void cartesian_debug(Point point)
define _val_constant(a, k, l, m)((const double) _constant[a.i -_NVARMAX]) @define val_diagonal(a
scalar cartesian_init_vertex_scalar(scalar s, const char *name)
static scalar * list_add_depends(scalar *list, scalar s)
#define boundary(...)
void init_const_vector(vector v, const char *name, double *val)
int bid
tensor cartesian_init_tensor(tensor t, const char *name)
bid new_bid()
void(* boundary_face)(vectorl)
scalar * list_clone(scalar *l)
static void write_stencil_index(int *index)
This displays a (1D,2D,3D) stencil index.
define double double char flags
void stencil_val_a(Point p, scalar s, int i, int j, int k, bool input, const char *file, int line)
vector cartesian_init_vector(vector v, const char *name)
static double periodic_bc(Point point, Point neighbor, scalar s, bool *data)
define double double _z
scalar new_const_scalar(const char *name, int i, double val)
BoundaryFunc default_vector_bc[]
void stencil_val(Point p, scalar s, int i, int j, int k, const char *file, int line, bool overflow)
void cartesian_boundary_level(scalar *list, int l)
static void cartesian_scalar_clone(scalar clone, scalar src)
static void periodic_boundary(int d)
void(* boundary_level)(scalar *, int l)
tensor init_symmetric_tensor(tensor t, const char *name)
void cartesian_boundary_face(vectorl vl)
static double interpolate_linear(Point point, scalar v, double xp=0., double yp=0., double zp=0.)
scalar cartesian_init_scalar(scalar s, const char *name)
void cartesian_methods()
static double antisymmetry(Point point, Point neighbor, scalar s, bool *data)
static void debug_plot(FILE *fp, const char *name, const char *cells, const char *stencil)
static double symmetry(Point point, Point neighbor, scalar s, bool *data)
define l
vector cartesian_init_face_vector(vector v, const char *name)
static char * replace_(const char *vname)
define double double char Reduce reductions
void output_cells(FILE *fp=stdout, coord c={0}, double size=0.)
void boundary_flux(vector *list) __attribute__((deprecated))
define m((k)==0 &&(l)==0 &&(m)==0) macro2 foreach_point(double _x=0.
trace void interpolate_array(scalar *list, coord *a, int n, double *v, bool linear=false)
Point locate(double xp=0, double yp=0, double zp=0)
Definition cartesian.h:358
define neighbor(o, p, q)((Point)
Definition cartesian.h:146
double real
Definition cartesian.h:6
if TRASH undef trash define trash(list) reset(list
free(list1)
define VT _attribute[s.i] v y scalar * list
Definition cartesian.h:276
int y
Definition common.h:76
#define define
int x
Definition common.h:76
int z
Definition common.h:76
vector(* init_vector)(vector, const char *)
Definition common.h:349
scalar(* init_vertex_scalar)(scalar, const char *)
Definition common.h:348
void(* scalar_clone)(scalar, scalar)
Definition common.h:352
scalar * all
Definition common.h:342
const scalar cm[]
Definition common.h:397
#define nodata
Definition common.h:7
static _Attributes * _attribute
Definition common.h:165
scalar * list_add(scalar *list, scalar s)
Definition common.h:208
tensor(* init_tensor)(tensor, const char *)
Definition common.h:351
void(* free_solver_func)(void)
Definition common.h:508
static Array * free_solver_funcs
Definition common.h:510
struct @0 Period
scalar * baseblock
Definition common.h:343
scalar * list_copy(scalar *l)
Definition common.h:225
void(* BoundaryStencilFunc)(Point, Point, scalar, void *)
Definition common.h:149
vector(* init_face_vector)(vector, const char *)
Definition common.h:350
static int sign(number x)
Definition common.h:13
scalar * list_append(scalar *list, scalar s)
Definition common.h:188
static bool is_vertex_scalar(scalar s)
Definition common.h:334
scalar(* init_scalar)(scalar, const char *)
Definition common.h:347
@ bottom
Definition common.h:123
@ left
Definition common.h:123
int nboundary
Definition common.h:127
double(* BoundaryFunc)(Point, Point, scalar, bool *)
Definition common.h:148
double * _constant
Definition common.h:131
size_t datasize
Definition common.h:132
scalar f[]
The primary fields are:
Definition two-phase.h:56
vector w[]
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line line line type
Definition config.h:571
static void qpclose_all()
Definition config.h:654
#define p
Definition tree.h:320
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line calloc(n, s) @ define prealloc(p
define m m double _val_higher_dimension
Definition config.h:587
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define file
Definition config.h:120
#define assert(a)
Definition config.h:107
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line realloc(p, s) @ define pfree(p
define sysmalloc malloc define syscalloc calloc define sysrealloc realloc define sysfree free define systrdup strdup define line line line line strdup(s) @ define tracing(...) @ define end_tracing(...) @define tid() 0 @define pid() 0 @define npe() 1 @define mpi_all_reduce(v
Point point
Definition conserving.h:86
coord o
Definition curvature.h:672
else return n
Definition curvature.h:101
trace bool box(bool notics=false, float lc[3]={0}, float lw=1.)
Definition draw.h:1397
trace bool cells(coord n={0, 0, 1}, double alpha=0., float lc[3]={0}, float lw=1.)
Definition draw.h:1128
scalar s
Definition embed-tree.h:56
macro2 double neumann(double expr, Point point=point, scalar s=_s, bool *data=data)
Definition embed.h:732
*cs[i, 0, 0] a *[i -1, 0, 0] j
Definition embed.h:88
macro2 double neumann_homogeneous(double expr, Point point=point, scalar s=_s, bool *data=data)
Definition embed.h:740
double v[2]
Definition embed.h:383
double d[2]
Definition embed.h:383
scalar int i
Definition embed.h:74
macro2 double dirichlet(double expr, Point point=point, scalar s=_s, bool *data=data)
For ease of use, we replace the Neumann and Dirichlet functions with macros so that they can be used ...
Definition embed.h:717
macro2 double dirichlet_homogeneous(double expr, Point point=point, scalar s=_s, bool *data=data)
Definition embed.h:725
static int line
Definition errors.c:754
static char * fname
Definition errors.c:753
double t
Definition events.h:24
static Event * Events
Definition events.h:21
#define realloc_scalar(size)
Definition gpu.h:480
#define free_grid()
Definition grid.h:1404
static int input(void)
Definition include.c:2085
static int
Definition include.c:978
#define kh_init(name)
Definition khash.h:430
khint_t khiter_t
Definition khash.h:154
#define kh_get(name, h, k)
Definition khash.h:474
#define kh_destroy(name, h)
Definition khash.h:437
#define khash_t(name)
Definition khash.h:423
#define kh_value(h, x)
Definition khash.h:512
#define KHASH_MAP_INIT_INT64(name, khval_t)
Definition khash.h:599
#define kh_end(h)
Definition khash.h:526
#define kh_put(name, h, k, r)
Definition khash.h:465
#define data(k, l)
Definition linear.h:26
void pmuntrace(void)
Definition mtrace.h:192
FILE * fp
Definition mtrace.h:7
size_t size
size *double * b
int int int level
char dir[]
Definition qcc.c:64
void boundary_stencil(ForeachData *loop)
This functions applies the boundary conditions, as defined by check_stencil().
Definition stencils.h:274
def _i
Definition stencils.h:405
def _k
Definition stencils.h:405
@ s_input
Definition stencils.h:23
@ s_centered
Definition stencils.h:19
@ s_output
Definition stencils.h:24
@ s_nowarning
Definition stencils.h:25
@ s_face
Definition stencils.h:20
static bool scalar_is_dirty(scalar s)
Definition stencils.h:74
def _j
Definition stencils.h:405
void * p
Definition array.h:6
long len
Definition array.h:7
Definition linear.h:21
int level
Definition linear.h:23
char * name
Definition common.h:156
Definition events.h:6
Event * next
Definition events.h:18
int last
Definition events.h:7
char * name
Definition stencils.h:42
Definition common.h:78
double x
Definition common.h:79
int i
Definition common.h:44
vector x
Definition common.h:67
scalar x
Definition common.h:47
scalar * x
Definition common.h:57
define n sizeof(Cell))) @define fine(a
#define periodic(dir)
Definition tree.h:1734
def allocated(k, l, n)(mem_allocated(((Tree *) grid) -> L[point.level]->m, point.i+k, point.j+l)) @ @def NEIGHBOR(k, l, n)(((((Tree *) grid) ->L[point.level]->m) ->b[point.i+k][point.j+l])) @ @def PARENT(k, l, n)(((((Tree *) grid) ->L[point.level-1]->m) ->b[(point.i+2)/2+k][(point.j+2)/2+l])) @ @def allocated_child(k, l, n)(level< depth() &&mem_allocated(((Tree *) grid) ->L[point.level+1]->m, 2 *point.i- 2+k, 2 *point.j- 2+l)) @ @def CHILD(k, l, n)(((((Tree *) grid) ->L[point.level+1]->m) ->b[2 *point.i- 2+k][2 *point.j- 2+l])) @ @define CELL(m)(*((Cell *)(m))) @define depth()(grid->depth) @define aparent(k, l, n) CELL(PARENT(k, l, n)) @define child(k, l, n) CELL(CHILD(k, l, n)) @define cell CELL(NEIGHBOR(0, 0, 0)) @define neighbor(k, l, n) CELL(NEIGHBOR(k, l, n)) @def neighborp(l, m, n)(Point)
Definition tree.h:253
Array * index
void * pointer
Definition view.h:321
#define lambda
scalar flux[]
Definition vof.h:166
scalar c
Definition vof.h:57
src vof fflush(ferr)