Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
rpe.h
Go to the documentation of this file.
1/** @file rpe.h
2 */
3/**
4# Resting Potential Energy
5
6See [Ilicak et al., 2012](#ilicak2012), appendix A and [Petersen et
7al., 2015](#petersen2015), section 2.3. This should obviously be much
8better documented than this.
9
10Note that other more computationally-efficient options exist, based
11for example on sampling of histograms of density.
12
13~~~bib
14@article{ilicak2012,
15title = {Spurious dianeutral mixing and the role of momentum closure},
16journal = {Ocean Modelling},
17volume = {45-46},
18pages = {37-58},
19year = {2012},
20issn = {1463-5003},
21doi = {https://doi.org/10.1016/j.ocemod.2011.10.003},
22url = {https://www.sciencedirect.com/science/article/pii/S1463500311001685},
23author = {Mehmet Ilicak and Alistair J. Adcroft and Stephen M. Griffies and
24Robert W. Hallberg}
25}
26
27@article{petersen2015,
28title = {Evaluation of the arbitrary {L}agrangian–{E}ulerian vertical coordinate
29 method in the {MPAS}-Ocean model},
30journal = {Ocean Modelling},
31volume = {86},
32pages = {93-113},
33year = {2015},
34issn = {1463-5003},
35doi = {https://doi.org/10.1016/j.ocemod.2014.12.004},
36url = {https://www.sciencedirect.com/science/article/pii/S1463500314001796},
37author = {Mark R. Petersen and Douglas W. Jacobsen and Todd D. Ringler and
38Matthew W. Hecht and Mathew E. Maltrud},
39}
40~~~
41*/
42
43#include "utils.h"
44#if dimension == 2
45# include "fractions.h"
46#endif
47
48struct _Zarea {
49 scalar zb; // compulsory
50 double * A, * V; // compulsory
51 int n; // compulsory
52 double max; // optional (default 0.)
53 double min; // optional (default zb.min)
54};
55typedef struct _Zarea Zarea;
56
57#define A(i) (a->A[i])
58
59double area_integral (Zarea * a, double z1, double z2, double * Az)
60{
61 assert (z2 >= z1);
62 double dz = (a->max - a->min)/(a->n - 1);
63
64 int i2 = (z2 - a->min)/(a->max - a->min)*(a->n - 1);
65 double f1, f2, z;
66 if (i2 < 0)
67 z = z2 = 0.;
68 else if (i2 < a->n - 1) {
69 f1 = A(i2), z = a->min + dz*i2;
70 f2 = A(i2) + (z2 - z)/dz*(A(i2 + 1) - f1);
71 }
72 else {
73 z = a->max, f1 = f2 = A(a->n - 1);
74 i2 = a->n - 1;
75 }
76 double v2, zv2;
77 if (z < z2) {
78 v2 = (z2 - z)*(f2 + f1)/2.;
79 zv2 = ((sq(z2) - sq(z))/2.*(f1*z2 - f2*z) +
80 (cube(z2) - cube(z))/3.*(f2 - f1))/(z2 - z);
81 }
82 else
83 v2 = 0., zv2 = 0.;
84
85 int i1;
86 double v1, zv1;
87 if (z1 <= a->min)
88 v1 = 0., zv1 = 0., i1 = -1;
89 else {
90 i1 = (z1 - a->min)/(a->max - a->min)*(a->n - 1);
91 double f1, f2, z;
92 if (i1 >= a->n - 1) {
93 f1 = f2 = A(a->n - 1);
94 z = a->max;
95 }
96 else {
97 f1 = A(i1) + (z1 - a->min - i1*dz)/dz*(A(i1 + 1) - A(i1));
98 if (i2 > i1)
99 f2 = A(i1 + 1), z = a->min + dz*(i1 + 1);
100 else
101 f2 = A(i1), z = a->min + dz*i1;
102 }
103 if (z != z1) {
104 v1 = (z - z1)*(f2 + f1)/2.;
105 zv1 = ((sq(z) - sq(z1))/2.*(f1*z - f2*z1) +
106 (cube(z) - cube(z1))/3.*(f2 - f1))/(z - z1);
107 }
108 else
109 v1 = 0., zv1 = 0.;
110 }
111
112 double V = v1 + v2;
113 *Az = zv1 + zv2;
114 for (int i = i1 + 1; i < i2; i++) {
115 z = a->min + i*dz, z2 = z + dz;
116 f1 = A(i), f2 = A(i+1);
117 V += (z2 - z)*(f2 + f1)/2.;
118 *Az += ((sq(z2) - sq(z))/2.*(f1*z2 - f2*z) +
119 (cube(z2) - cube(z))/3.*(f2 - f1))/(z2 - z);
120 }
121 return V;
122}
123
124double area_z2 (Zarea * a, double z1, double dv, double * Az)
125{
126 double zm = z1, z2 = z1 + 10.*(a->max - a->min);
127 double vol = area_integral (a, z1, z2, Az);
128 assert (vol > dv);
129 while (fabs (z2 - zm) > dry) {
130 double z = (zm + z2)/2.;
131 vol = area_integral (a, z1, z, Az);
132 if (vol > dv) z2 = z;
133 else zm = z;
134 // fprintf (stderr, "^^^ dv = %g: zm: %g z2: %g vol: %g\n", dv, zm, z2, vol);
135 }
136 // assert (fabs(dv - vol) < 1e-6*dv);
137 return (zm + z2)/2.;
138}
139
140Zarea zarea (scalar zb, double * A, double * V, int n,
141 double max = 0., double min = 0.)
142{
143 if (min == 0.)
144 min = statsf (zb).min;
145 for (int i = 0; i < n; i++) {
146 double val = min + i*(max - min)/(n - 1);
147#if dimension == 2
148 scalar f[];
149 fractions (zb, f, val = val);
150 stats sf = statsf(f);
151 A[i] = sf.volume - sf.sum;
152#else
153 double area = 0.;
154 foreach(reduction(+:area)) {
155 double a = (zb[] + zb[-1])/2. - val, b = (zb[] + zb[1])/2. - val;
156 if (a*b < 0.) {
157 double x = a/(a - b);
158 area += dv()*(a < 0. ? x : 1. - x);
159 }
160 else if (a < 0.)
161 area += dv();
162 }
163 A[i] = area;
164#endif
165 }
166 return (Zarea){ zb, A, V, n, max, min };
167}
168
169Zarea zvolume (scalar zb, double * A, double * V, int n,
170 double max = 0., double min = 0.)
171{
172 Zarea s = zarea (zb, A, V, n, max, min);
173 double volume = 0.;
174 for (int i = 0; i < n; i++) {
175 double dz = (s.max - s.min)/(n - 1);
176 volume += dz*A[i];
177 V[i] = volume;
178 }
179 return s;
180}
181
182double volumez (Zarea s, double volume)
183{
184 if (volume <= 0.) return s.min;
185 if (volume >= s.V[s.n-1]) return s.max;
186 int a = 0, b = s.n - 1;
187 while (b > a + 1) {
188 int i = (a + b)/2;
189 if (s.V[i] > volume) b = i;
190 else a = i;
191 }
192 double c = (s.V[a] - volume)/(s.V[a] - s.V[b]);
193 double dz = (s.max - s.min)/(s.n - 1);
194 return (s.min + a*dz)*(1. - c) + (s.min + b*dz)*c;
195}
196
197double zareaval (Zarea s, double z)
198{
199 double dz = (s.max - s.min)/(s.n - 1);
200 int j = (z - s.min)/dz;
201 assert (j >= 0 && j < s.n - 1);
202#if 0
203 printf ("A %d %g %g %g\n", j, s.A[j], s.A[j+1],
204 s.A[j] + (z - s.min)/dz*(s.A[j+1] - s.A[j]));
205#endif
206 return s.A[j] + (z - s.min)/dz*(s.A[j+1] - s.A[j]);
207}
208
209double barycenter (Zarea s, double z1, double z2)
210{
211 if (z1 == z2)
212 return z1;
213 assert (z1 >= s.min && z1 < z2 && z2 <= s.max);
214 int n = 1;
215 double dz = (z2 - z1)/n, sz = 0., sa = 0., z = z1 + dz/2.;
216 for (int i = 0; i < n; i++) {
217 double a = zareaval (s, z);
218 sa += a, sz += z*a, z += dz;
219 }
220 // printf ("# %g %g %g\n", z1, z2, sz/sa);
221 return sz/sa;
222}
223
224typedef struct {
225 double rho, dv, z;
226} Parcel;
227
228int heavier_than (const void * a, const void * b)
229{
230 Parcel * p1 = (Parcel *)a, * p2 = (Parcel *)b;
231 return p1->rho > p2->rho ? -1 : 1;
232}
233
234double energy (double * PE = NULL, double * KE = NULL)
235{
236 double vPE = 0., vKE = 0.;
237 foreach(reduction(+:vPE) reduction(+:vKE)) {
238 double z = zb[];
239 foreach_layer() {
240 z += h[]/2;
241 double mass = h[]*dv()*(1. + drho(T[]));
242 vPE += mass*z;
243#if NH
244 vKE += mass*sq(w[])/2.;
245#endif
246 for (int _d = 0; _d < dimension; _d++)
247 vKE += mass*sq(u.x[])/2.;
248 z += h[]/2;
249 }
250 }
251 vPE *= G;
252 if (PE) *PE = vPE;
253 if (KE) *KE = vKE;
254 return vPE + vKE;
255}
256
257// Parallel sort
258
259size_t psort (void ** base, size_t nmemb, size_t size,
260 int (* compar)(const void *, const void *))
261{
262#if _MPI // fixme: sorting is not done in parallel
263 // Number of MPI processes and current rank
264 int nproc, rank;
267
268 int counts[nproc];
269 // Each process tells the root how many elements it holds
270 nmemb *= size;
272
273 // Displacements in the receive buffer for MPI_GATHERV
274 int disps[nproc];
275 // Displacement for the first chunk of data - 0
276 for (int i = 0; i < nproc; i++)
277 disps[i] = i > 0 ? (disps[i - 1] + counts[i - 1]) : 0;
278
279 if (rank == 0) {
280 size_t nmemb1 = disps[nproc-1] + counts[nproc-1];
281 void * alldata = malloc (nmemb1);
284 free (*base);
285 *base = alldata;
286 nmemb = nmemb1/size;
287 qsort (*base, nmemb, size, compar);
288 }
289 else
292#else
293 qsort (*base, nmemb, size, compar);
294#endif
295 return nmemb;
296}
297
298// Resting Potential Energy: see e.g. Ilicak et al, 2012, Appendix A
299
300trace
301double RPE (int n = 100)
302{
303 double A[n], V[n];
304 Zarea vol = zarea (zb, A, V, n);
305#if 0
306 for (int i = 0; i < vol.n; i++) {
307 double dz = (vol.max - vol.min)/(vol.n - 1);
308 fprintf (stderr, "%g %g %g\n", vol.min + i*dz, V[i], - statsf (zb).sum);
309 }
310
311 for (double volume = 0.; volume <= V[vol.n-1]; volume += V[vol.n-1]/35.)
312 printf ("%g %g\n", volumez (vol, volume), volume);
313 exit (0);
314#endif
315
316 long nb = 0;
317 for (int _i = 0; _i < _N; _i++) /* foreach_leaf */ nb++;
318 nb *= nl;
319 Parcel * p = malloc (sizeof(Parcel)*nb);
320 long i = 0;
321 for (int _i = 0; _i < _N; _i++) /* foreach_leaf */
322 foreach_layer() {
323 assert (i < nb);
324 p[i++] = (Parcel){ 1. + drho(T[]), h[]*dv() };
325 }
326
327 nb = psort ((void **) &p, nb, sizeof (Parcel), heavier_than);
328
329 double RPE = HUGE;
330 if (pid() == 0) {
331 RPE = 0.;
332 double volume = 0.;
333 double z = vol.min;
334 for (int i = 0; i < nb; i++)
335 if (p[i].dv > 0.) {
336 volume += p[i].dv;
337#if 1
338 double Az, zt = area_z2 (&vol, z, p[i].dv, &Az);
339 p[i].z = Az/p[i].dv;
340 RPE += p[i].rho*Az; // see e.g. (1) in Ilicak et al, 2012
341#else
342 double zt = volumez (vol, volume);
343 p[i].z = barycenter (vol, z, zt);
344
345 // fprintf (stderr, "&& %g %g %g %g\n", zt, zt1, p[i].z, Az/p[i].dv);
346
347 p[i].z = Az/p[i].dv;
348 RPE += p[i].rho*p[i].z*p[i].dv; // see e.g. (1) in Ilicak et al, 2012
349#endif
350
351 // fprintf (stdout, "%g %g %g %g\n", p[i].rho, p[i].dv, p[i].z, zt);
352 z = zt;
353 }
354 }
355 // fprintf (stdout, "****\n");
357 free (p);
358
359 return RPE*G;
360}
#define u
Definition advection.h:30
const vector a
Definition all-mach.h:59
scalar zb[]
Definition atmosphere.h:6
double G
Definition atmosphere.h:12
scalar h[]
Definition atmosphere.h:6
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
free(list1)
int x
Definition common.h:76
int z
Definition common.h:76
#define HUGE
Definition common.h:6
static number sq(number x)
Definition common.h:11
#define dv()
Definition common.h:92
static number cube(number x)
Definition common.h:12
scalar f[]
The primary fields are:
Definition two-phase.h:56
vector w[]
#define p
Definition tree.h:320
#define assert(a)
Definition config.h:107
else return n
Definition curvature.h:101
scalar s
Definition embed-tree.h:56
*cs[i, 0, 0] a *[i -1, 0, 0] j
Definition embed.h:88
scalar int i
Definition embed.h:74
double dry
Definition entrainment.h:30
trace void fractions(scalar Phi, scalar c, vector s={0}, double val=0.)
Definition fractions.h:123
double * drho
static float area(const KdtRect rect)
Definition kdt.c:439
double * sum
int nl
Definition layers.h:9
#define foreach_layer()
#define rho(f)
The density and viscosity are defined using arithmetic averages by default.
Definition momentum.h:45
size_t max
Definition mtrace.h:8
size_t size
size *double * b
int heavier_than(const void *a, const void *b)
Definition rpe.h:228
size_t psort(void **base, size_t nmemb, size_t size, int(*compar)(const void *, const void *))
Definition rpe.h:259
Zarea zvolume(scalar zb, double *A, double *V, int n, double max=0., double min=0.)
Definition rpe.h:169
double barycenter(Zarea s, double z1, double z2)
Definition rpe.h:209
double area_z2(Zarea *a, double z1, double dv, double *Az)
Definition rpe.h:124
Zarea zarea(scalar zb, double *A, double *V, int n, double max=0., double min=0.)
Definition rpe.h:140
#define A(i)
Definition rpe.h:57
double volumez(Zarea s, double volume)
Definition rpe.h:182
double zareaval(Zarea s, double z)
Definition rpe.h:197
double area_integral(Zarea *a, double z1, double z2, double *Az)
Definition rpe.h:59
trace double RPE(int n=100)
Definition rpe.h:301
double energy(double *PE=NULL, double *KE=NULL)
Definition rpe.h:234
def _i
Definition stencils.h:405
Definition rpe.h:224
double dv
Definition rpe.h:225
Definition rpe.h:48
double min
Definition rpe.h:53
double max
Definition rpe.h:52
scalar zb
Definition rpe.h:49
int n
Definition rpe.h:51
double * A
Definition rpe.h:50
double * V
Definition rpe.h:50
The statsf() function returns the minimum, maximum, volume sum, standard deviation and volume for fie...
Definition utils.h:166
double min
Definition utils.h:167
scalar T[]
Definition thermal.h:60
#define sf
We have the option of using some "smearing" of the density/viscosity jump.
stats statsf(scalar f)
Definition utils.h:170
scalar c
Definition vof.h:57