-
Notifications
You must be signed in to change notification settings - Fork 6
Expand file tree
/
Copy pathcalc_drho.cpp
More file actions
477 lines (336 loc) · 12.5 KB
/
calc_drho.cpp
File metadata and controls
477 lines (336 loc) · 12.5 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
#include <petscksp.h>
#include <petscdmda.h>
extern Vec dRho;
extern Vec Temper;
extern Vec geoq_rho;
extern double alpha_exp_thermo;
extern double gravity;
extern double RHOM;
PetscErrorCode calc_drho()
{
PetscErrorCode ierr=0;
VecCopy(geoq_rho,dRho);
VecScale(dRho,-gravity);
// VecCopy(Temper, dRho); CHKERRQ(ierr);
// PetscScalar a;
// a=alpha_exp_thermo*gravity*RHOM; CHKERRQ(ierr);//check: replace RHOM by geoq_rho
// VecAXPBY(dRho,-gravity,a,geoq_rho); CHKERRQ(ierr);
return ierr;
}
extern Vec local_P;
extern Vec local_P_aux;
extern Vec Pressure;
extern Vec Pressure_aux;
extern Vec local_geoq_rho;
extern DM da_Veloc;
extern DM da_Thermal;
extern long Nx,Ny,Nz;
extern long V_GT;
extern long T_NE;
extern double Lx, Ly, depth;
extern int n_interfaces;
extern PetscScalar *interfaces;
extern PetscScalar *inter_rho;
typedef struct {
PetscScalar u;
PetscScalar w;
} Stokes2d;
typedef struct {
PetscScalar u;
PetscScalar v;
PetscScalar w;
} Stokes3d;
PetscErrorCode DMDAGetElementCorners_2d(DM da,PetscInt *sx,PetscInt *sz,PetscInt *mx,PetscInt *mz);
PetscErrorCode DMDAGetElementCorners_3d(DM da,PetscInt *sx,PetscInt *sy, PetscInt *sz,PetscInt *mx, PetscInt *my, PetscInt *mz);
PetscErrorCode calc_pressure_2d()
{
Stokes2d **pp;
PetscErrorCode ierr;
ierr = VecZeroEntries(local_P);CHKERRQ(ierr);
ierr = DMGlobalToLocalBegin(da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMGlobalToLocalEnd( da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMDAVecGetArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
int M,P;
double zz;
PetscFunctionBeginUser;
ierr = DMDAGetInfo(da_Thermal,0,&M,&P,NULL,0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
PetscInt sex1,sez1,mx1,mz1;
PetscInt sex,sez,mx,mz;
PetscInt ei,ek;
ierr = DMDAGetElementCorners_2d(da_Thermal,&sex,&sez,&mx,&mz);CHKERRQ(ierr);
ierr = DMDAGetElementCorners_2d(da_Veloc,&sex1,&sez1,&mx1,&mz1);CHKERRQ(ierr);
if ((sex-sex1!=0)||(sez-sez1!=0)||(mx-mx1!=0)||(mz-mz1!=0)){
printf("%d %d %d %d\n",sex-sex1,sez-sez1,mx-mx1,mz-mz1);
SETERRQ1(PETSC_COMM_WORLD,1,"Wrong partition (temper,velocity)\n",1);
}
PetscReal interp_interfaces[n_interfaces];
for (ek = sez; ek < sez+mz; ek++) {
for (ei = sex; ei < sex+mx; ei++) {
/*indr[0].i=ei ; indr[0].j=ek ;
indr[1].i=ei+1; indr[1].j=ek ;
indr[2].i=ei ; indr[2].j=ek+1;
indr[3].i=ei+1; indr[3].j=ek+1;
xx = ei*Lx/(M-1);*/
zz = -(P-1-ek)*depth/(P-1);
for (int in=0;in<n_interfaces;in++){
float rfac = (0.5);
interp_interfaces[in] = interfaces[ei + Nx*in] * rfac;
rfac = (0.5);
interp_interfaces[in] += interfaces[(ei+1) + Nx*in] * rfac;
}
double pressure_cumulat = 0;
int verif=0;
/*for (int in=0;in<n_interfaces && verif==0;in++){
if (zz<interp_interfaces[in]){
verif=1;
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-zz);
}
else{
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in+1]-interp_interfaces[in]);
}
}
if (verif==0){
pressure_cumulat += inter_rho[n_interfaces]*gravity*(-zz);
}*/
/*int in;
for (in=n_interfaces; in>0 && verif==0;in--){
if (zz>=interp_interfaces[in]){
verif=1;
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-zz);
}
else{
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in+1]-interp_interfaces[in]);
}
}
if (verif==0){
pressure_cumulat += inter_rho[0]*gravity*(-zz);
}*/
if (zz>interp_interfaces[n_interfaces-1]){
pressure_cumulat += inter_rho[n_interfaces]*gravity*(-zz);
}
else{
pressure_cumulat += inter_rho[n_interfaces]*gravity*
(-interp_interfaces[n_interfaces-1]);
for (int in=n_interfaces-1; in>0 && verif==0;in--){
if (zz>interp_interfaces[in-1]){
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-zz);
verif=1;
}
else{
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-interp_interfaces[in-1]);
}
}
if (verif==0){
pressure_cumulat += inter_rho[0]*gravity*(interp_interfaces[0]-zz);
}
}
pp[ek][ei].u = pressure_cumulat;
}
}
ierr = DMDAVecRestoreArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
return ierr;
}
PetscErrorCode shift_pressure_2d() //necessary if the surface pressure is not close to zero
{
Stokes2d **pp;
PetscScalar **pp_aux;
PetscErrorCode ierr;
//get Pressure array
ierr = DMGlobalToLocalBegin(da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMGlobalToLocalEnd( da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMDAVecGetArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
//get Pressure_aux array
ierr = DMGlobalToLocalBegin(da_Thermal,Pressure_aux,INSERT_VALUES,local_P_aux);
ierr = DMGlobalToLocalEnd( da_Thermal,Pressure_aux,INSERT_VALUES,local_P_aux);
ierr = DMDAVecGetArray(da_Thermal,local_P_aux,&pp_aux);CHKERRQ(ierr);
PetscInt sx,sz,mmx,mmz;
PetscInt i,k;
ierr = DMDAGetCorners(da_Thermal,&sx,&sz,NULL,&mmx,&mmz,NULL);CHKERRQ(ierr);
PetscScalar ppp;
PetscInt cont_ppp;
for (k=sz; k<sz+mmz; k++) {
for (i=sx; i<sx+mmx; i++) {
ppp=0.0;
cont_ppp = 0;
if (k<Nz-1 && i<Nx-1) {ppp+=pp[k ][i ].u; cont_ppp++;}
if (k<Nz-1 && i>0 ) {ppp+=pp[k ][i-1].u; cont_ppp++;}
if (k>0 && i<Nx-1) {ppp+=pp[k-1][i ].u; cont_ppp++;}
if (k>0 && i>0 ) {ppp+=pp[k-1][i-1].u; cont_ppp++;}
pp_aux[k][i] = ppp/cont_ppp;
}
}
int cont_air = 0, cont_air_all;
float pressure_air = 0.0, pressure_air_all;
for (k=sz; k<sz+mmz; k++) {
for (i=sx; i<sx+mmx; i++) {
if (k==Nz-1){
pressure_air+=pp_aux[k][i];
cont_air++;
}
}
}
//restore Pressure_aux
ierr = DMDAVecRestoreArray(da_Thermal,local_P_aux,&pp_aux);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Thermal,local_P_aux,INSERT_VALUES,Pressure_aux);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Thermal,local_P_aux,INSERT_VALUES,Pressure_aux);CHKERRQ(ierr);
MPI_Allreduce(&cont_air,&cont_air_all,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
MPI_Allreduce(&pressure_air,&pressure_air_all,1,MPI_FLOAT,MPI_SUM,PETSC_COMM_WORLD);
PetscReal pressure_min;
pressure_min = pressure_air_all/cont_air_all;
pressure_min=-pressure_min;
VecShift(Pressure_aux,pressure_min);
//restore Pressure
ierr = DMDAVecRestoreArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
return ierr;
}
PetscErrorCode calc_pressure_3d()
{
Stokes3d ***pp;
PetscErrorCode ierr;
ierr = VecZeroEntries(local_P);CHKERRQ(ierr);
ierr = DMGlobalToLocalBegin(da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMGlobalToLocalEnd( da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMDAVecGetArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
//PetscScalar **qq_rho;
//ierr = DMGlobalToLocalBegin(da_Thermal,geoq_rho,INSERT_VALUES,local_geoq_rho);
//ierr = DMGlobalToLocalEnd( da_Thermal,geoq_rho,INSERT_VALUES,local_geoq_rho);
//ierr = DMDAVecGetArray(da_Thermal,local_geoq_rho,&qq_rho);CHKERRQ(ierr);
int M,N,P;
double zz;
PetscFunctionBeginUser;
ierr = DMDAGetInfo(da_Thermal,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);CHKERRQ(ierr);
printf("%d %d %d P\n",M,N,P);
//
PetscInt sex1,sey1,sez1,mx1,my1,mz1;
PetscInt sex,sey,sez,mx,my,mz;
PetscInt ei,ej,ek;
ierr = DMDAGetElementCorners_3d(da_Thermal,&sex,&sey,&sez,&mx,&my,&mz);CHKERRQ(ierr);
ierr = DMDAGetElementCorners_3d(da_Veloc,&sex1,&sey1,&sez1,&mx1,&my1,&mz1);CHKERRQ(ierr);
if ((sex-sex1!=0)||(sey-sey1!=0)||(sez-sez1!=0)||(mx-mx1!=0)||(my-my1!=0)||(mz-mz1!=0)){
printf("%d %d %d %d %d %d\n",sex-sex1,sey-sey1,sez-sez1,mx-mx1,my-my1,mz-mz1);
SETERRQ1(PETSC_COMM_WORLD,1,"Wrong partition (temper,velocity)\n",1);
}
PetscReal interp_interfaces[n_interfaces];
for (ek = sez; ek < sez+mz; ek++) {
for (ej = sey; ej < sey+my; ej++) {
for (ei = sex; ei < sex+mx; ei++) {
zz = -(P-1-ek)*depth/(P-1);
for (int in=0;in<n_interfaces;in++){
float rfac = 0.25;
interp_interfaces[in] = interfaces[ej*Nx+ei + Nx*Ny*in] * rfac;
rfac = 0.25;
interp_interfaces[in] += interfaces[ej*Nx+(ei+1) + Nx*Ny*in] * rfac;
rfac = 0.25;
interp_interfaces[in] += interfaces[(ej+1)*Nx+ei + Nx*Ny*in] * rfac;
rfac = 0.25;
interp_interfaces[in] += interfaces[(ej+1)*Nx+(ei+1) + Nx*Ny*in] * rfac;
}
double pressure_cumulat = 0;
int verif=0;
if (zz>interp_interfaces[n_interfaces-1]){
pressure_cumulat += inter_rho[n_interfaces]*gravity*(-zz);
}
else{
pressure_cumulat += inter_rho[n_interfaces]*gravity*
(-interp_interfaces[n_interfaces-1]);
for (int in=n_interfaces-1; in>0 && verif==0;in--){
if (zz>interp_interfaces[in-1]){
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-zz);
verif=1;
}
else{
pressure_cumulat += inter_rho[in]*gravity*
(interp_interfaces[in]-interp_interfaces[in-1]);
}
}
if (verif==0){
pressure_cumulat += inter_rho[0]*gravity*(interp_interfaces[0]-zz);
}
}
pp[ek][ej][ei].u = pressure_cumulat;
//printf("%d ",in);
//if (ek==1) printf("pp = %lf\n",pressure_cumulat);
}
}
}
//ierr = DMDAVecRestoreArray(da_Thermal,local_geoq_rho,&qq_rho);CHKERRQ(ierr);
ierr = DMDAVecRestoreArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
//write_pressure(-1);
return ierr;
}
PetscErrorCode shift_pressure_3d() //necessary if the surface pressure is not close to zero
{
PetscErrorCode ierr;
Stokes3d ***pp;
PetscScalar ***pp_aux;
//get Pressure array
ierr = DMGlobalToLocalBegin(da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMGlobalToLocalEnd( da_Veloc,Pressure,INSERT_VALUES,local_P);
ierr = DMDAVecGetArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
//get Pressure_aux array
ierr = DMGlobalToLocalBegin(da_Thermal,Pressure_aux,INSERT_VALUES,local_P_aux);
ierr = DMGlobalToLocalEnd( da_Thermal,Pressure_aux,INSERT_VALUES,local_P_aux);
ierr = DMDAVecGetArray(da_Thermal,local_P_aux,&pp_aux);CHKERRQ(ierr);
PetscInt sx,sy,sz,mmx,mmy,mmz;
PetscInt i,j,k;
ierr = DMDAGetCorners(da_Veloc,&sx,&sy,&sz,&mmx,&mmy,&mmz);CHKERRQ(ierr);
PetscScalar ppp;
PetscInt cont_ppp;
for (k=sz; k<sz+mmz; k++) {
for (j=sy; j<sy+mmy; j++) {
for (i=sx; i<sx+mmx; i++) {
ppp=0.0;
cont_ppp = 0;
if (k<Nz-1 && j<Ny-1 && i<Nx-1) {ppp+=pp[k ][j ][i ].u; cont_ppp++;}
if (k<Nz-1 && j<Ny-1 && i>0 ) {ppp+=pp[k ][j ][i-1].u; cont_ppp++;}
if (k<Nz-1 && j>0 && i<Nx-1) {ppp+=pp[k ][j-1][i ].u; cont_ppp++;}
if (k<Nz-1 && j>0 && i>0 ) {ppp+=pp[k ][j-1][i-1].u; cont_ppp++;}
if (k>0 && j<Ny-1 && i<Nx-1) {ppp+=pp[k-1][j ][i ].u; cont_ppp++;}
if (k>0 && j<Ny-1 && i>0 ) {ppp+=pp[k-1][j ][i-1].u; cont_ppp++;}
if (k>0 && j>0 && i<Nx-1) {ppp+=pp[k-1][j-1][i ].u; cont_ppp++;}
if (k>0 && j>0 && i>0 ) {ppp+=pp[k-1][j-1][i-1].u; cont_ppp++;}
pp_aux[k][j][i] = ppp/cont_ppp;
}
}
}
int cont_air = 0, cont_air_all;
float pressure_air = 0.0, pressure_air_all;
for (k=sz; k<sz+mmz; k++) {
for (j=sy; j<sy+mmy; j++) {
for (i=sx; i<sx+mmx; i++) {
if (k==Nz-1){
pressure_air+=pp_aux[k][j][i];
cont_air++;
}
}
}
}
//restore Pressure_aux
ierr = DMDAVecRestoreArray(da_Thermal,local_P_aux,&pp_aux);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Thermal,local_P_aux,INSERT_VALUES,Pressure_aux);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Thermal,local_P_aux,INSERT_VALUES,Pressure_aux);CHKERRQ(ierr);
MPI_Allreduce(&cont_air,&cont_air_all,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);
MPI_Allreduce(&pressure_air,&pressure_air_all,1,MPI_FLOAT,MPI_SUM,PETSC_COMM_WORLD);
PetscReal pressure_min;
pressure_min = pressure_air_all/cont_air_all;
pressure_min=-pressure_min;
VecShift(Pressure_aux,pressure_min);
//restore Pressure
ierr = DMDAVecRestoreArray(da_Veloc,local_P,&pp);CHKERRQ(ierr);
ierr = DMLocalToGlobalBegin(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
ierr = DMLocalToGlobalEnd(da_Veloc,local_P,INSERT_VALUES,Pressure);CHKERRQ(ierr);
return ierr;
}