1 #ifndef AMREX_EB2_3D_C_H_
2 #define AMREX_EB2_3D_C_H_
3 #include <AMReX_Config.H>
22 [=] (
int i,
int j,
int k) noexcept
24 if ( s(i,j ,k ) < 0.0_rt && s(i+1,j ,k ) < 0.0_rt
25 && s(i,j+1,k ) < 0.0_rt && s(i+1,j+1,k ) < 0.0_rt
26 && s(i,j ,k+1) < 0.0_rt && s(i+1,j ,k+1) < 0.0_rt
27 && s(i,j+1,k+1) < 0.0_rt && s(i+1,j+1,k+1) < 0.0_rt)
29 cell(i,j,k).setRegular();
31 else if (s(i,j ,k ) >= 0.0_rt && s(i+1,j ,k ) >= 0.0_rt
32 && s(i,j+1,k ) >= 0.0_rt && s(i+1,j+1,k ) >= 0.0_rt
33 && s(i,j ,k+1) >= 0.0_rt && s(i+1,j ,k+1) >= 0.0_rt
34 && s(i,j+1,k+1) >= 0.0_rt && s(i+1,j+1,k+1) >= 0.0_rt)
36 cell(i,j,k).setCovered();
40 cell(i,j,k).setSingleValued();
49 [=] (
int i,
int j,
int k) noexcept
51 if ( s(i,j,k ) < 0.0_rt && s(i,j+1,k ) < 0.0_rt
52 && s(i,j,k+1) < 0.0_rt && s(i,j+1,k+1) < 0.0_rt )
56 else if (s(i,j,k ) >= 0.0_rt && s(i,j+1,k ) >= 0.0_rt
57 && s(i,j,k+1) >= 0.0_rt && s(i,j+1,k+1) >= 0.0_rt )
72 [=] (
int i,
int j,
int k) noexcept
74 if ( s(i,j,k ) < 0.0_rt && s(i+1,j,k ) < 0.0_rt
75 && s(i,j,k+1) < 0.0_rt && s(i+1,j,k+1) < 0.0_rt )
79 else if (s(i,j,k ) >= 0.0_rt && s(i+1,j,k ) >= 0.0_rt
80 && s(i,j,k+1) >= 0.0_rt && s(i+1,j,k+1) >= 0.0_rt )
95 [=] (
int i,
int j,
int k) noexcept
97 if ( s(i,j ,k) < 0.0_rt && s(i+1,j ,k) < 0.0_rt
98 && s(i,j+1,k) < 0.0_rt && s(i+1,j+1,k) < 0.0_rt)
102 else if (s(i,j ,k) >= 0.0_rt && s(i+1,j ,k) >= 0.0_rt
103 && s(i,j+1,k) >= 0.0_rt && s(i+1,j+1,k) >= 0.0_rt)
118 [=] (
int i,
int j,
int k) noexcept
120 if (s(i,j,k) < 0.0_rt && s(i+1,j,k) < 0.0_rt) {
122 }
else if (s(i,j,k) >= 0.0_rt && s(i+1,j,k) >= 0.0_rt) {
134 [=] (
int i,
int j,
int k) noexcept
136 if (s(i,j,k) < 0.0_rt && s(i,j+1,k) < 0.0_rt) {
138 }
else if (s(i,j,k) >= 0.0_rt && s(i,j+1,k) >= 0.0_rt) {
150 [=] (
int i,
int j,
int k) noexcept
152 if (s(i,j,k) < 0.0_rt && s(i,j,k+1) < 0.0_rt) {
154 }
else if (s(i,j,k) >= 0.0_rt && s(i,j,k+1) >= 0.0_rt) {
164 int num_cuts (Real a, Real b) noexcept {
165 return (a >= 0.0_rt && b < 0.0_rt) || (
b >= 0.0_rt && a < 0.0_rt);
170 int check_mvmc (
int i,
int j,
int k, Array4<Real const>
const&
fine)
181 int nx00 =
num_cuts(
fine(i,j,k),
fine(i+1,j,k)) +
num_cuts(
fine(i+1,j,k),
fine(i+2,j,k));
182 int nx10 =
num_cuts(
fine(i,j+2,k),
fine(i+1,j+2,k)) +
num_cuts(
fine(i+1,j+2,k),
fine(i+2,j+2,k));
183 int nx01 =
num_cuts(
fine(i,j,k+2),
fine(i+1,j,k+2)) +
num_cuts(
fine(i+1,j,k+2),
fine(i+2,j,k+2));
184 int nx11 =
num_cuts(
fine(i,j+2,k+2),
fine(i+1,j+2,k+2)) +
num_cuts(
fine(i+1,j+2,k+2),
fine(i+2,j+2,k+2));
187 int ny00 =
num_cuts(
fine(i,j,k),
fine(i,j+1,k)) +
num_cuts(
fine(i,j+1,k),
fine(i,j+2,k));
188 int ny10 =
num_cuts(
fine(i+2,j,k),
fine(i+2,j+1,k)) +
num_cuts(
fine(i+2,j+1,k),
fine(i+2,j+2,k));
189 int ny01 =
num_cuts(
fine(i,j,k+2),
fine(i,j+1,k+2)) +
num_cuts(
fine(i,j+1,k+2),
fine(i,j+2,k+2));
190 int ny11 =
num_cuts(
fine(i+2,j,k+2),
fine(i+2,j+1,k+2)) +
num_cuts(
fine(i+2,j+1,k+2),
fine(i+2,j+2,k+2));
193 int nz00 =
num_cuts(
fine(i,j,k),
fine(i,j,k+1)) +
num_cuts(
fine(i,j,k+1),
fine(i,j,k+2));
194 int nz10 =
num_cuts(
fine(i+2,j,k),
fine(i+2,j,k+1)) +
num_cuts(
fine(i+2,j,k+1),
fine(i+2,j,k+2));
195 int nz01 =
num_cuts(
fine(i,j+2,k),
fine(i,j+2,k+1)) +
num_cuts(
fine(i,j+2,k+1),
fine(i,j+2,k+2));
196 int nz11 =
num_cuts(
fine(i+2,j+2,k),
fine(i+2,j+2,k+1)) +
num_cuts(
fine(i+2,j+2,k+1),
fine(i+2,j+2,k+2));
200 int n = ny00 + ny01 + nz00 + nz01;
210 n = ny10 + ny11 + nz10 + nz11;
221 n = nx00 + nx01 + nz00 + nz10;
231 n = nx10 + nx11 + nz01 + nz11;
242 n = nx00 + nx10 + ny00 + ny10;
252 n = nx01 + nx11 + ny01 + ny11;
261 if (nxm == 1 && nym == 1 && nzm == 1 && nxp == 1 && nyp == 1 && nzp == 1) {
262 n = (
fine(i ,j ,k ) < 0.0_rt) + (
fine(i+2,j ,k ) < 0.0_rt) +
263 (
fine(i ,j+2,k ) < 0.0_rt) + (
fine(i+2,j+2,k ) < 0.0_rt) +
264 (
fine(i ,j ,k+2) < 0.0_rt) + (
fine(i+2,j ,k+2) < 0.0_rt) +
265 (
fine(i ,j+2,k+2) < 0.0_rt) + (
fine(i+2,j+2,k+2) < 0.0_rt);
266 if (n == 2 || n == 6) {
270 amrex::Abort(
"amrex::check_mvmc: how did this happen? nopen != 4");
281 if (f1 == Real(1.0) && f2 == Real(1.0)) {
283 }
else if (f1 == Real(-1.0) && f2 == Real(-1.0)) {
286 if (f1 == Real(-1.0)) {
288 }
else if (f1 == Real(1.0)) {
291 f1 = Real(0.5)*f1 - Real(0.25);
293 if (f2 == Real(-1.0)) {
295 }
else if (f2 == Real(1.0)) {
298 f2 = Real(0.5)*f2 + Real(0.25);
300 Real
r = (f2*f2-f1*f1)/(f2-f1+Real(1.e-30));
348 if (fflag(ii,jj ,kk ).isRegular() && fflag(ii+1,jj ,kk ).isRegular() &&
349 fflag(ii,jj+1,kk ).isRegular() && fflag(ii+1,jj+1,kk ).isRegular() &&
350 fflag(ii,jj ,kk+1).isRegular() && fflag(ii+1,jj ,kk+1).isRegular() &&
351 fflag(ii,jj+1,kk+1).isRegular() && fflag(ii+1,jj+1,kk+1).isRegular())
353 cflag(i,j,k).setRegular();
354 cvol(i,j,k) = 1.0_rt;
355 ccent(i,j,k,0) = 0.0_rt;
356 ccent(i,j,k,1) = 0.0_rt;
357 ccent(i,j,k,2) = 0.0_rt;
359 cbc(i,j,k,0) = -1.0_rt;
360 cbc(i,j,k,1) = -1.0_rt;
361 cbc(i,j,k,2) = -1.0_rt;
362 cbn(i,j,k,0) = 0.0_rt;
363 cbn(i,j,k,1) = 0.0_rt;
364 cbn(i,j,k,2) = 0.0_rt;
366 else if (fflag(ii,jj ,kk ).isCovered() && fflag(ii+1,jj ,kk ).isCovered() &&
367 fflag(ii,jj+1,kk ).isCovered() && fflag(ii+1,jj+1,kk ).isCovered() &&
368 fflag(ii,jj ,kk+1).isCovered() && fflag(ii+1,jj ,kk+1).isCovered() &&
369 fflag(ii,jj+1,kk+1).isCovered() && fflag(ii+1,jj+1,kk+1).isCovered())
371 cflag(i,j,k).setCovered();
372 cvol(i,j,k) = 0.0_rt;
373 ccent(i,j,k,0) = 0.0_rt;
374 ccent(i,j,k,1) = 0.0_rt;
375 ccent(i,j,k,2) = 0.0_rt;
377 cbc(i,j,k,0) = -1.0_rt;
378 cbc(i,j,k,1) = -1.0_rt;
379 cbc(i,j,k,2) = -1.0_rt;
380 cbn(i,j,k,0) = 0.0_rt;
381 cbn(i,j,k,1) = 0.0_rt;
382 cbn(i,j,k,2) = 0.0_rt;
386 cflag(i,j,k).setSingleValued();
388 cvol(i,j,k) = 0.125_rt*(fvol(ii ,jj ,kk ) + fvol(ii+1,jj ,kk ) +
389 fvol(ii ,jj+1,kk ) + fvol(ii+1,jj+1,kk ) +
390 fvol(ii ,jj ,kk+1) + fvol(ii+1,jj ,kk+1) +
391 fvol(ii ,jj+1,kk+1) + fvol(ii+1,jj+1,kk+1));
392 Real cvolinv = 1.0_rt / cvol(i,j,k);
394 ccent(i,j,k,0) = 0.125_rt * cvolinv *
395 (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,0)-0.25_rt) +
396 fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,0)+0.25_rt) +
397 fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,0)-0.25_rt) +
398 fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,0)+0.25_rt) +
399 fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,0)-0.25_rt) +
400 fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,0)+0.25_rt) +
401 fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,0)-0.25_rt) +
402 fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,0)+0.25_rt));
403 ccent(i,j,k,1) = 0.125_rt * cvolinv *
404 (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,1)-0.25_rt) +
405 fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,1)-0.25_rt) +
406 fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,1)+0.25_rt) +
407 fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,1)+0.25_rt) +
408 fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,1)-0.25_rt) +
409 fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,1)-0.25_rt) +
410 fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,1)+0.25_rt) +
411 fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,1)+0.25_rt));
412 ccent(i,j,k,2) = 0.125_rt * cvolinv *
413 (fvol(ii ,jj ,kk )*(0.5_rt*fcent(ii ,jj ,kk ,2)-0.25_rt) +
414 fvol(ii+1,jj ,kk )*(0.5_rt*fcent(ii+1,jj ,kk ,2)-0.25_rt) +
415 fvol(ii ,jj+1,kk )*(0.5_rt*fcent(ii ,jj+1,kk ,2)-0.25_rt) +
416 fvol(ii+1,jj+1,kk )*(0.5_rt*fcent(ii+1,jj+1,kk ,2)-0.25_rt) +
417 fvol(ii ,jj ,kk+1)*(0.5_rt*fcent(ii ,jj ,kk+1,2)+0.25_rt) +
418 fvol(ii+1,jj ,kk+1)*(0.5_rt*fcent(ii+1,jj ,kk+1,2)+0.25_rt) +
419 fvol(ii ,jj+1,kk+1)*(0.5_rt*fcent(ii ,jj+1,kk+1,2)+0.25_rt) +
420 fvol(ii+1,jj+1,kk+1)*(0.5_rt*fcent(ii+1,jj+1,kk+1,2)+0.25_rt));
422 cba(i,j,k) = 0.25_rt*(fba(ii ,jj ,kk ) + fba(ii+1,jj ,kk ) +
423 fba(ii ,jj+1,kk ) + fba(ii+1,jj+1,kk ) +
424 fba(ii ,jj ,kk+1) + fba(ii+1,jj ,kk+1) +
425 fba(ii ,jj+1,kk+1) + fba(ii+1,jj+1,kk+1));
426 Real cbainv = 1.0_rt / cba(i,j,k);
428 cbc(i,j,k,0) = 0.25_rt * cbainv *
429 ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,0)-0.25_rt)
430 + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,0)+0.25_rt)
431 + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,0)-0.25_rt)
432 + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,0)+0.25_rt)
433 + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,0)-0.25_rt)
434 + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,0)+0.25_rt)
435 + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,0)-0.25_rt)
436 + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,0)+0.25_rt) );
437 cbc(i,j,k,1) = 0.25_rt * cbainv *
438 ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,1)-0.25_rt)
439 + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,1)-0.25_rt)
440 + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,1)+0.25_rt)
441 + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,1)+0.25_rt)
442 + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,1)-0.25_rt)
443 + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,1)-0.25_rt)
444 + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,1)+0.25_rt)
445 + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,1)+0.25_rt) );
446 cbc(i,j,k,2) = 0.25_rt * cbainv *
447 ( fba(ii ,jj ,kk )*(0.5_rt*fbc(ii ,jj ,kk ,2)-0.25_rt)
448 + fba(ii+1,jj ,kk )*(0.5_rt*fbc(ii+1,jj ,kk ,2)-0.25_rt)
449 + fba(ii ,jj+1,kk )*(0.5_rt*fbc(ii ,jj+1,kk ,2)-0.25_rt)
450 + fba(ii+1,jj+1,kk )*(0.5_rt*fbc(ii+1,jj+1,kk ,2)-0.25_rt)
451 + fba(ii ,jj ,kk+1)*(0.5_rt*fbc(ii ,jj ,kk+1,2)+0.25_rt)
452 + fba(ii+1,jj ,kk+1)*(0.5_rt*fbc(ii+1,jj ,kk+1,2)+0.25_rt)
453 + fba(ii ,jj+1,kk+1)*(0.5_rt*fbc(ii ,jj+1,kk+1,2)+0.25_rt)
454 + fba(ii+1,jj+1,kk+1)*(0.5_rt*fbc(ii+1,jj+1,kk+1,2)+0.25_rt) );
456 Real nx = fbn(ii ,jj ,kk ,0)*fba(ii ,jj ,kk )
457 + fbn(ii+1,jj ,kk ,0)*fba(ii+1,jj ,kk )
458 + fbn(ii ,jj+1,kk ,0)*fba(ii ,jj+1,kk )
459 + fbn(ii+1,jj+1,kk ,0)*fba(ii+1,jj+1,kk )
460 + fbn(ii ,jj ,kk+1,0)*fba(ii ,jj ,kk+1)
461 + fbn(ii+1,jj ,kk+1,0)*fba(ii+1,jj ,kk+1)
462 + fbn(ii ,jj+1,kk+1,0)*fba(ii ,jj+1,kk+1)
463 + fbn(ii+1,jj+1,kk+1,0)*fba(ii+1,jj+1,kk+1);
464 Real ny = fbn(ii ,jj ,kk ,1)*fba(ii ,jj ,kk )
465 + fbn(ii+1,jj ,kk ,1)*fba(ii+1,jj ,kk )
466 + fbn(ii ,jj+1,kk ,1)*fba(ii ,jj+1,kk )
467 + fbn(ii+1,jj+1,kk ,1)*fba(ii+1,jj+1,kk )
468 + fbn(ii ,jj ,kk+1,1)*fba(ii ,jj ,kk+1)
469 + fbn(ii+1,jj ,kk+1,1)*fba(ii+1,jj ,kk+1)
470 + fbn(ii ,jj+1,kk+1,1)*fba(ii ,jj+1,kk+1)
471 + fbn(ii+1,jj+1,kk+1,1)*fba(ii+1,jj+1,kk+1);
472 Real nz = fbn(ii ,jj ,kk ,2)*fba(ii ,jj ,kk )
473 + fbn(ii+1,jj ,kk ,2)*fba(ii+1,jj ,kk )
474 + fbn(ii ,jj+1,kk ,2)*fba(ii ,jj+1,kk )
475 + fbn(ii+1,jj+1,kk ,2)*fba(ii+1,jj+1,kk )
476 + fbn(ii ,jj ,kk+1,2)*fba(ii ,jj ,kk+1)
477 + fbn(ii+1,jj ,kk+1,2)*fba(ii+1,jj ,kk+1)
478 + fbn(ii ,jj+1,kk+1,2)*fba(ii ,jj+1,kk+1)
479 + fbn(ii+1,jj+1,kk+1,2)*fba(ii+1,jj+1,kk+1);
480 Real nfac = 1.0_rt /
std::sqrt(nx*nx+ny*ny+nz*nz+1.e-30_rt);
481 cbn(i,j,k,0) = nx*nfac;
482 cbn(i,j,k,1) = ny*nfac;
483 cbn(i,j,k,2) = nz*nfac;
484 ierr = (nx == 0.0_rt && ny == 0.0_rt && nz == 0.0_rt)
487 || ( fapx(ii,jj ,kk )==1.0_rt && fapx(ii+2,jj ,kk )==1.0_rt
488 && fapx(ii,jj+1,kk )==1.0_rt && fapx(ii+2,jj+1,kk )==1.0_rt
489 && fapx(ii,jj ,kk+1)==1.0_rt && fapx(ii+2,jj ,kk+1)==1.0_rt
490 && fapx(ii,jj+1,kk+1)==1.0_rt && fapx(ii+2,jj+1,kk+1)==1.0_rt
491 && fapy(ii,jj ,kk )==1.0_rt && fapy(ii+1,jj ,kk )==1.0_rt
492 && fapy(ii,jj+2,kk )==1.0_rt && fapy(ii+1,jj+2,kk )==1.0_rt
493 && fapy(ii,jj ,kk+1)==1.0_rt && fapy(ii+1,jj ,kk+1)==1.0_rt
494 && fapy(ii,jj+2,kk+1)==1.0_rt && fapy(ii+1,jj+2,kk+1)==1.0_rt
495 && fapz(ii,jj ,kk )==1.0_rt && fapz(ii+1,jj ,kk )==1.0_rt
496 && fapz(ii,jj+1,kk )==1.0_rt && fapz(ii+1,jj+1,kk )==1.0_rt
497 && fapz(ii,jj ,kk+2)==1.0_rt && fapz(ii+1,jj ,kk+2)==1.0_rt
498 && fapz(ii,jj+1,kk+2)==1.0_rt && fapz(ii+1,jj+1,kk+2)==1.0_rt);
504 cvol(i,j,k) = 1.0_rt;
505 ccent(i,j,k,0) = 0.0_rt;
506 ccent(i,j,k,1) = 0.0_rt;
507 ccent(i,j,k,2) = 0.0_rt;
509 cbc(i,j,k,0) = -1.0_rt;
510 cbc(i,j,k,1) = -1.0_rt;
511 cbc(i,j,k,2) = -1.0_rt;
512 cbn(i,j,k,0) = 0.0_rt;
513 cbn(i,j,k,1) = 0.0_rt;
514 cbn(i,j,k,2) = 0.0_rt;
519 capx(i,j,k) = 0.25_rt*(fapx(ii,jj ,kk ) + fapx(ii,jj+1,kk ) +
520 fapx(ii,jj ,kk+1) + fapx(ii,jj+1,kk+1));
521 if (capx(i,j,k) != 0.0_rt) {
522 Real apinv = 1.0_rt / capx(i,j,k);
523 cfcx(i,j,k,0) = 0.25_rt * apinv *
524 (fapx(ii,jj ,kk )*(0.5_rt*ffcx(ii,jj ,kk ,0)-0.25_rt) +
525 fapx(ii,jj+1,kk )*(0.5_rt*ffcx(ii,jj+1,kk ,0)+0.25_rt) +
526 fapx(ii,jj ,kk+1)*(0.5_rt*ffcx(ii,jj ,kk+1,0)-0.25_rt) +
527 fapx(ii,jj+1,kk+1)*(0.5_rt*ffcx(ii,jj+1,kk+1,0)+0.25_rt) );
528 cfcx(i,j,k,1) = 0.25_rt * apinv *
529 (fapx(ii,jj ,kk )*(0.5_rt*ffcx(ii,jj ,kk ,1)-0.25_rt) +
530 fapx(ii,jj+1,kk )*(0.5_rt*ffcx(ii,jj+1,kk ,1)-0.25_rt) +
531 fapx(ii,jj ,kk+1)*(0.5_rt*ffcx(ii,jj ,kk+1,1)+0.25_rt) +
532 fapx(ii,jj+1,kk+1)*(0.5_rt*ffcx(ii,jj+1,kk+1,1)+0.25_rt) );
534 cfcx(i,j,k,0) = 0.0_rt;
535 cfcx(i,j,k,1) = 0.0_rt;
540 capx(i,j,k) = 1.0_rt;
541 cfcx(i,j,k,0) = 0.0_rt;
542 cfcx(i,j,k,1) = 0.0_rt;
547 capy(i,j,k) = 0.25_rt*(fapy(ii ,jj,kk ) + fapy(ii+1,jj,kk ) +
548 fapy(ii ,jj,kk+1) + fapy(ii+1,jj,kk+1));
549 if (capy(i,j,k) != 0.0_rt) {
550 Real apinv = 1.0_rt / capy(i,j,k);
551 cfcy(i,j,k,0) = 0.25_rt * apinv *
552 (fapy(ii ,jj,kk )*(0.5_rt*ffcy(ii ,jj,kk ,0)-0.25_rt) +
553 fapy(ii+1,jj,kk )*(0.5_rt*ffcy(ii+1,jj,kk ,0)+0.25_rt) +
554 fapy(ii ,jj,kk+1)*(0.5_rt*ffcy(ii ,jj,kk+1,0)-0.25_rt) +
555 fapy(ii+1,jj,kk+1)*(0.5_rt*ffcy(ii+1,jj,kk+1,0)+0.25_rt) );
556 cfcy(i,j,k,1) = 0.25_rt * apinv *
557 (fapy(ii ,jj,kk )*(0.5_rt*ffcy(ii ,jj,kk ,1)-0.25_rt) +
558 fapy(ii+1,jj,kk )*(0.5_rt*ffcy(ii+1,jj,kk ,1)-0.25_rt) +
559 fapy(ii ,jj,kk+1)*(0.5_rt*ffcy(ii ,jj,kk+1,1)+0.25_rt) +
560 fapy(ii+1,jj,kk+1)*(0.5_rt*ffcy(ii+1,jj,kk+1,1)+0.25_rt) );
562 cfcy(i,j,k,0) = 0.0_rt;
563 cfcy(i,j,k,1) = 0.0_rt;
568 capy(i,j,k) = 1.0_rt;
569 cfcy(i,j,k,0) = 0.0_rt;
570 cfcy(i,j,k,1) = 0.0_rt;
575 capz(i,j,k) = 0.25_rt * (fapz(ii ,jj ,kk) + fapz(ii+1,jj ,kk) +
576 fapz(ii ,jj+1,kk) + fapz(ii+1,jj+1,kk));
577 if (capz(i,j,k) != 0.0_rt) {
578 Real apinv = 1.0_rt / capz(i,j,k);
579 cfcz(i,j,k,0) = 0.25_rt * apinv *
580 (fapz(ii ,jj ,kk)*(0.5_rt*ffcz(ii ,jj ,kk,0)-0.25_rt) +
581 fapz(ii+1,jj ,kk)*(0.5_rt*ffcz(ii+1,jj ,kk,0)+0.25_rt) +
582 fapz(ii ,jj+1,kk)*(0.5_rt*ffcz(ii ,jj+1,kk,0)-0.25_rt) +
583 fapz(ii+1,jj+1,kk)*(0.5_rt*ffcz(ii+1,jj+1,kk,0)+0.25_rt) );
584 cfcz(i,j,k,1) = 0.25_rt * apinv *
585 (fapz(ii ,jj ,kk)*(0.5_rt*ffcz(ii ,jj ,kk,1)-0.25_rt) +
586 fapz(ii+1,jj ,kk)*(0.5_rt*ffcz(ii+1,jj ,kk,1)-0.25_rt) +
587 fapz(ii ,jj+1,kk)*(0.5_rt*ffcz(ii ,jj+1,kk,1)+0.25_rt) +
588 fapz(ii+1,jj+1,kk)*(0.5_rt*ffcz(ii+1,jj+1,kk,1)+0.25_rt) );
590 cfcz(i,j,k,0) = 0.0_rt;
591 cfcz(i,j,k,1) = 0.0_rt;
596 capz(i,j,k) = 1.0_rt;
597 cfcz(i,j,k,0) = 0.0_rt;
598 cfcz(i,j,k,1) = 0.0_rt;
607 cecx(i,j,k) = 1.0_rt;
616 cecy(i,j,k) = 1.0_rt;
625 cecz(i,j,k) = 1.0_rt;
636 auto flg = cflag(i,j,k);
637 flg.setDisconnected();
639 if (!flg.isCovered())
641 flg.setConnected(0,0,0);
643 if (apx(i ,j,k) != 0.0_rt) { flg.setConnected(-1, 0, 0); }
644 if (apx(i+1,j,k) != 0.0_rt) { flg.setConnected( 1, 0, 0); }
645 if (apy(i,j ,k) != 0.0_rt) { flg.setConnected( 0, -1, 0); }
646 if (apy(i,j+1,k) != 0.0_rt) { flg.setConnected( 0, 1, 0); }
647 if (apz(i,j,k ) != 0.0_rt) { flg.setConnected( 0, 0, -1); }
648 if (apz(i,j,k+1) != 0.0_rt) { flg.setConnected( 0, 0, 1); }
650 if ( (apx(i,j,k) != 0.0_rt && apy(i-1,j,k) != 0.0_rt) ||
651 (apy(i,j,k) != 0.0_rt && apx(i,j-1,k) != 0.0_rt) )
653 flg.setConnected(-1, -1, 0);
654 if (apz(i-1,j-1,k ) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
655 if (apz(i-1,j-1,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
658 if ( (apx(i+1,j,k) != 0.0_rt && apy(i+1,j ,k) != 0.0_rt) ||
659 (apy(i ,j,k) != 0.0_rt && apx(i+1,j-1,k) != 0.0_rt) )
661 flg.setConnected(1, -1, 0);
662 if (apz(i+1,j-1,k ) != 0.0_rt) { flg.setConnected(1,-1,-1); }
663 if (apz(i+1,j-1,k+1) != 0.0_rt) { flg.setConnected(1,-1, 1); }
666 if ( (apx(i,j ,k) != 0.0_rt && apy(i-1,j+1,k) != 0.0_rt) ||
667 (apy(i,j+1,k) != 0.0_rt && apx(i ,j+1,k) != 0.0_rt) )
669 flg.setConnected(-1, 1, 0);
670 if (apz(i-1,j+1,k ) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
671 if (apz(i-1,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
674 if ( (apx(i+1,j ,k) != 0.0_rt && apy(i+1,j+1,k) != 0.0_rt) ||
675 (apy(i ,j+1,k) != 0.0_rt && apx(i+1,j+1,k) != 0.0_rt) )
677 flg.setConnected(1, 1, 0);
678 if (apz(i+1,j+1,k ) != 0.0_rt) { flg.setConnected(1, 1,-1); }
679 if (apz(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected(1, 1, 1); }
682 if ( (apx(i,j,k) != 0.0_rt && apz(i-1,j,k ) != 0.0_rt) ||
683 (apz(i,j,k) != 0.0_rt && apx(i ,j,k-1) != 0.0_rt) )
685 flg.setConnected(-1, 0, -1);
686 if (apy(i-1,j ,k-1) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
687 if (apy(i-1,j+1,k-1) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
690 if ( (apx(i+1,j,k) != 0.0_rt && apz(i+1,j,k ) != 0.0_rt) ||
691 (apz(i ,j,k) != 0.0_rt && apx(i+1,j,k-1) != 0.0_rt) )
693 flg.setConnected(1, 0, -1);
694 if (apy(i+1,j ,k-1) != 0.0_rt) { flg.setConnected(1,-1,-1); }
695 if (apy(i+1,j+1,k-1) != 0.0_rt) { flg.setConnected(1, 1,-1); }
698 if ( (apx(i,j,k ) != 0.0_rt && apz(i-1,j,k+1) != 0.0_rt) ||
699 (apz(i,j,k+1) != 0.0_rt && apx(i ,j,k+1) != 0.0_rt) )
701 flg.setConnected(-1, 0, 1);
702 if (apy(i-1,j ,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
703 if (apy(i-1,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
706 if ( (apx(i+1,j,k ) != 0.0_rt && apz(i+1,j,k+1) != 0.0_rt) ||
707 (apz(i ,j,k+1) != 0.0_rt && apx(i+1,j,k+1) != 0.0_rt) )
709 flg.setConnected(1, 0, 1);
710 if (apy(i+1,j ,k+1) != 0.0_rt) { flg.setConnected(1,-1, 1); }
711 if (apy(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected(1, 1, 1); }
714 if ( (apy(i,j,k) != 0.0_rt && apz(i,j-1,k ) != 0.0_rt) ||
715 (apz(i,j,k) != 0.0_rt && apy(i,j ,k-1) != 0.0_rt) )
717 flg.setConnected(0, -1, -1);
718 if (apx(i ,j-1,k-1) != 0.0_rt) { flg.setConnected(-1,-1,-1); }
719 if (apx(i+1,j-1,k-1) != 0.0_rt) { flg.setConnected( 1,-1,-1); }
722 if ( (apy(i,j+1,k) != 0.0_rt && apz(i,j+1,k ) != 0.0_rt) ||
723 (apz(i,j ,k) != 0.0_rt && apy(i,j+1,k-1) != 0.0_rt) )
725 flg.setConnected(0, 1, -1);
726 if (apx(i ,j+1,k-1) != 0.0_rt) { flg.setConnected(-1, 1,-1); }
727 if (apx(i+1,j+1,k-1) != 0.0_rt) { flg.setConnected( 1, 1,-1); }
730 if ( (apy(i,j,k ) != 0.0_rt && apz(i,j-1,k+1) != 0.0_rt) ||
731 (apz(i,j,k+1) != 0.0_rt && apy(i,j ,k+1) != 0.0_rt) )
733 flg.setConnected(0, -1, 1);
734 if (apx(i ,j-1,k+1) != 0.0_rt) { flg.setConnected(-1,-1, 1); }
735 if (apx(i+1,j-1,k+1) != 0.0_rt) { flg.setConnected( 1,-1, 1); }
738 if ( (apy(i,j+1,k ) != 0.0_rt && apz(i,j+1,k+1) != 0.0_rt) ||
739 (apz(i,j ,k+1) != 0.0_rt && apy(i,j+1,k+1) != 0.0_rt) )
741 flg.setConnected(0, 1, 1);
742 if (apx(i ,j+1,k+1) != 0.0_rt) { flg.setConnected(-1, 1, 1); }
743 if (apx(i+1,j+1,k+1) != 0.0_rt) { flg.setConnected( 1, 1, 1); }
#define AMREX_FORCE_INLINE
Definition: AMReX_Extension.H:119
#define AMREX_GPU_HOST_DEVICE
Definition: AMReX_GpuQualifiers.H:20
Array4< Real > fine
Definition: AMReX_InterpFaceRegister.cpp:90
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition: AMReX_Box.H:204
static constexpr Type_t regular
Definition: AMReX_EB2_Graph.H:38
static constexpr Type_t covered
Definition: AMReX_EB2_Graph.H:39
static constexpr Type_t irregular
Definition: AMReX_EB2_Graph.H:40
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int num_cuts(Real a, Real b) noexcept
Definition: AMReX_EB2_2D_C.H:71
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real coarsen_edge_cent(Real f1, Real f2)
Definition: AMReX_EB2_3D_C.H:279
Definition: AMReX_FabArrayBase.H:32
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int coarsen_from_fine(int i, int j, Box const &bx, int ngrow, Array4< Real > const &cvol, Array4< Real > const &ccent, Array4< Real > const &cba, Array4< Real > const &cbc, Array4< Real > const &cbn, Array4< Real > const &capx, Array4< Real > const &capy, Array4< Real > const &cfcx, Array4< Real > const &cfcy, Array4< Real > const &cecx, Array4< Real > const &cecy, Array4< EBCellFlag > const &cflag, Array4< Real const > const &fvol, Array4< Real const > const &fcent, Array4< Real const > const &fba, Array4< Real const > const &fbc, Array4< Real const > const &fbn, Array4< Real const > const &fapx, Array4< Real const > const &fapy, Array4< Real const > const &ffcx, Array4< Real const > const &ffcy, Array4< Real const > const &fecx, Array4< Real const > const &fecy, Array4< EBCellFlag const > const &fflag)
Definition: AMReX_EB2_2D_C.H:96
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void build_cellflag_from_ap(int i, int j, Array4< EBCellFlag > const &cflag, Array4< Real const > const &apx, Array4< Real const > const &apy)
Definition: AMReX_EB2_2D_C.H:268
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void amrex_eb2_build_types(Box const &tbx, Box const &bxg2, Array4< Real const > const &s, Array4< EBCellFlag > const &cell, Array4< Type_t > const &fx, Array4< Type_t > const &fy)
Definition: AMReX_EB2_2D_C.H:9
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int check_mvmc(int i, int j, int, Array4< Real const > const &fine)
Definition: AMReX_EB2_2D_C.H:77
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 min_ubound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition: AMReX_Box.H:1928
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:35
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition: AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & min(const T &a, const T &b) noexcept
Definition: AMReX_Algorithm.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Returns a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition: AMReX_Box.H:1399
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition: AMReX.cpp:225
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Dim3 max_lbound(BoxND< dim > const &b1, BoxND< dim > const &b2) noexcept
Definition: AMReX_Box.H:1909
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE GpuComplex< T > sqrt(const GpuComplex< T > &a_z) noexcept
Return the square root of a complex number.
Definition: AMReX_GpuComplex.H:373
AMREX_GPU_HOST_DEVICE AMREX_ATTRIBUTE_FLATTEN_FOR void Loop(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition: AMReX_Loop.H:125
Definition: AMReX_FabArrayCommI.H:896
Definition: AMReX_Array4.H:61