Calculates the net flux across the immersed boundary surface.
Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/poisson.h.
2473{
2474 PetscErrorCode ierr;
2475
2476
2478
2479
2481
2482
2483 DM da = user->
da, fda = user->
fda;
2484
2485 DMDALocalInfo info = user->
info;
2486
2487 PetscInt xs = info.xs, xe = info.xs + info.xm;
2488 PetscInt ys = info.ys, ye = info.ys + info.ym;
2489 PetscInt zs = info.zs, ze = info.zs + info.zm;
2490 PetscInt mx = info.mx, my = info.my, mz = info.mz;
2491
2492 PetscInt i, j, k,ibi;
2493 PetscInt lxs, lys, lzs, lxe, lye, lze;
2494
2495 lxs = xs; lxe = xe;
2496 lys = ys; lye = ye;
2497 lzs = zs; lze = ze;
2498
2499 if (xs==0) lxs = xs+1;
2500 if (ys==0) lys = ys+1;
2501 if (zs==0) lzs = zs+1;
2502
2503 if (xe==mx) lxe = xe-1;
2504 if (ye==my) lye = ye-1;
2505 if (ze==mz) lze = ze-1;
2506
2507 PetscReal epsilon=1.e-8;
2508 PetscReal ***nvert, ibmval=1.9999;
2509
2510 struct Components {
2511 PetscReal x;
2512 PetscReal y;
2513 PetscReal z;
2514 }***ucor, ***lucor,***csi, ***eta, ***zet;
2515
2516
2517 PetscInt xend=mx-2 ,yend=my-2,zend=mz-2;
2518
2522
2523 DMDAVecGetArray(fda, user->
Ucont, &ucor);
2524 DMDAVecGetArray(fda, user->
lCsi, &csi);
2525 DMDAVecGetArray(fda, user->
lEta, &eta);
2526 DMDAVecGetArray(fda, user->
lZet, &zet);
2527 DMDAVecGetArray(da, user->
lNvert, &nvert);
2528
2529 PetscReal libm_Flux, libm_area, libm_Flux_abs=0., ibm_Flux_abs;
2530 libm_Flux = 0;
2531 libm_area = 0;
2532
2534
2535
2536 PetscReal *lIB_Flux = NULL, *lIB_area = NULL, *IB_Flux = NULL, *IB_Area = NULL;
2537 if (NumberOfBodies > 1) {
2538
2539 lIB_Flux=(PetscReal *)calloc(NumberOfBodies,sizeof(PetscReal));
2540 lIB_area=(PetscReal *)calloc(NumberOfBodies,sizeof(PetscReal));
2541 IB_Flux=(PetscReal *)calloc(NumberOfBodies,sizeof(PetscReal));
2542 IB_Area=(PetscReal *)calloc(NumberOfBodies,sizeof(PetscReal));
2543
2544
2545 for (ibi=0; ibi<NumberOfBodies; ibi++) {
2546 lIB_Flux[ibi]=0.0;
2547 lIB_area[ibi]=0.0;
2548 IB_Flux[ibi]=0.0;
2549 IB_Area[ibi]=0.0;
2550 }
2551 }
2552
2553
2554
2555
2556
2557
2559
2560 for (k=lzs; k<lze; k++) {
2561 for (j=lys; j<lye; j++) {
2562 for (i=lxs; i<lxe; i++) {
2563 if (nvert[k][j][i] < 0.1) {
2564 if (nvert[k][j][i+1] > 0.1 && nvert[k][j][i+1] < ibmval && i < xend) {
2565
2566 if (fabs(ucor[k][j][i].x)>epsilon) {
2567 libm_Flux += ucor[k][j][i].x;
2568 if (flg==3)
2569 libm_Flux_abs += fabs(ucor[k][j][i].x)/sqrt(csi[k][j][i].x * csi[k][j][i].x +
2570 csi[k][j][i].y * csi[k][j][i].y +
2571 csi[k][j][i].z * csi[k][j][i].z);
2572 else
2573 libm_Flux_abs += fabs(ucor[k][j][i].x);
2574
2575 libm_area += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2576 csi[k][j][i].y * csi[k][j][i].y +
2577 csi[k][j][i].z * csi[k][j][i].z);
2578
2579 if (NumberOfBodies > 1) {
2580
2581 ibi=(int)((nvert[k][j][i+1]-1.0)*1001);
2582 lIB_Flux[ibi] += ucor[k][j][i].x;
2583 lIB_area[ibi] += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2584 csi[k][j][i].y * csi[k][j][i].y +
2585 csi[k][j][i].z * csi[k][j][i].z);
2586 }
2587 } else
2588 ucor[k][j][i].x=0.;
2589
2590 }
2591 if (nvert[k][j+1][i] > 0.1 && nvert[k][j+1][i] < ibmval && j < yend) {
2592
2593 if (fabs(ucor[k][j][i].y)>epsilon) {
2594 libm_Flux += ucor[k][j][i].y;
2595 if (flg==3)
2596 libm_Flux_abs += fabs(ucor[k][j][i].y)/sqrt(eta[k][j][i].x * eta[k][j][i].x +
2597 eta[k][j][i].y * eta[k][j][i].y +
2598 eta[k][j][i].z * eta[k][j][i].z);
2599 else
2600 libm_Flux_abs += fabs(ucor[k][j][i].y);
2601 libm_area += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2602 eta[k][j][i].y * eta[k][j][i].y +
2603 eta[k][j][i].z * eta[k][j][i].z);
2604 if (NumberOfBodies > 1) {
2605
2606 ibi=(int)((nvert[k][j+1][i]-1.0)*1001);
2607
2608 lIB_Flux[ibi] += ucor[k][j][i].y;
2609 lIB_area[ibi] += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2610 eta[k][j][i].y * eta[k][j][i].y +
2611 eta[k][j][i].z * eta[k][j][i].z);
2612 }
2613 } else
2614 ucor[k][j][i].y=0.;
2615 }
2616 if (nvert[k+1][j][i] > 0.1 && nvert[k+1][j][i] < ibmval && k < zend) {
2617
2618 if (fabs(ucor[k][j][i].z)>epsilon) {
2619 libm_Flux += ucor[k][j][i].z;
2620 if (flg==3)
2621 libm_Flux_abs += fabs(ucor[k][j][i].z)/sqrt(zet[k][j][i].x * zet[k][j][i].x +
2622 zet[k][j][i].y * zet[k][j][i].y +
2623 zet[k][j][i].z * zet[k][j][i].z);
2624 else
2625 libm_Flux_abs += fabs(ucor[k][j][i].z);
2626 libm_area += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2627 zet[k][j][i].y * zet[k][j][i].y +
2628 zet[k][j][i].z * zet[k][j][i].z);
2629
2630 if (NumberOfBodies > 1) {
2631
2632 ibi=(int)((nvert[k+1][j][i]-1.0)*1001);
2633 lIB_Flux[ibi] += ucor[k][j][i].z;
2634 lIB_area[ibi] += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2635 zet[k][j][i].y * zet[k][j][i].y +
2636 zet[k][j][i].z * zet[k][j][i].z);
2637 }
2638 }else
2639 ucor[k][j][i].z=0.;
2640 }
2641 }
2642
2643 if (nvert[k][j][i] > 0.1 && nvert[k][j][i] < ibmval) {
2644
2645 if (nvert[k][j][i+1] < 0.1 && i < xend) {
2646 if (fabs(ucor[k][j][i].x)>epsilon) {
2647 libm_Flux -= ucor[k][j][i].x;
2648 if (flg==3)
2649 libm_Flux_abs += fabs(ucor[k][j][i].x)/sqrt(csi[k][j][i].x * csi[k][j][i].x +
2650 csi[k][j][i].y * csi[k][j][i].y +
2651 csi[k][j][i].z * csi[k][j][i].z);
2652 else
2653 libm_Flux_abs += fabs(ucor[k][j][i].x);
2654 libm_area += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2655 csi[k][j][i].y * csi[k][j][i].y +
2656 csi[k][j][i].z * csi[k][j][i].z);
2657 if (NumberOfBodies > 1) {
2658
2659 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2660 lIB_Flux[ibi] -= ucor[k][j][i].x;
2661 lIB_area[ibi] += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2662 csi[k][j][i].y * csi[k][j][i].y +
2663 csi[k][j][i].z * csi[k][j][i].z);
2664 }
2665
2666 }else
2667 ucor[k][j][i].x=0.;
2668 }
2669 if (nvert[k][j+1][i] < 0.1 && j < yend) {
2670 if (fabs(ucor[k][j][i].y)>epsilon) {
2671 libm_Flux -= ucor[k][j][i].y;
2672 if (flg==3)
2673 libm_Flux_abs += fabs(ucor[k][j][i].y)/ sqrt(eta[k][j][i].x * eta[k][j][i].x +
2674 eta[k][j][i].y * eta[k][j][i].y +
2675 eta[k][j][i].z * eta[k][j][i].z);
2676 else
2677 libm_Flux_abs += fabs(ucor[k][j][i].y);
2678 libm_area += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2679 eta[k][j][i].y * eta[k][j][i].y +
2680 eta[k][j][i].z * eta[k][j][i].z);
2681 if (NumberOfBodies > 1) {
2682
2683 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2684 lIB_Flux[ibi] -= ucor[k][j][i].y;
2685 lIB_area[ibi] += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2686 eta[k][j][i].y * eta[k][j][i].y +
2687 eta[k][j][i].z * eta[k][j][i].z);
2688 }
2689 }else
2690 ucor[k][j][i].y=0.;
2691 }
2692 if (nvert[k+1][j][i] < 0.1 && k < zend) {
2693 if (fabs(ucor[k][j][i].z)>epsilon) {
2694 libm_Flux -= ucor[k][j][i].z;
2695 if (flg==3)
2696 libm_Flux_abs += fabs(ucor[k][j][i].z)/sqrt(zet[k][j][i].x * zet[k][j][i].x +
2697 zet[k][j][i].y * zet[k][j][i].y +
2698 zet[k][j][i].z * zet[k][j][i].z);
2699 else
2700 libm_Flux_abs += fabs(ucor[k][j][i].z);
2701 libm_area += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2702 zet[k][j][i].y * zet[k][j][i].y +
2703 zet[k][j][i].z * zet[k][j][i].z);
2704 if (NumberOfBodies > 1) {
2705
2706 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2707 lIB_Flux[ibi] -= ucor[k][j][i].z;
2708 lIB_area[ibi] += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2709 zet[k][j][i].y * zet[k][j][i].y +
2710 zet[k][j][i].z * zet[k][j][i].z);
2711 }
2712 }else
2713 ucor[k][j][i].z=0.;
2714 }
2715 }
2716
2717 }
2718 }
2719 }
2720
2721 ierr = MPI_Allreduce(&libm_Flux, ibm_Flux,1,MPI_DOUBLE,MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
2722 ierr = MPI_Allreduce(&libm_Flux_abs, &ibm_Flux_abs,1,MPI_DOUBLE,MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
2723 ierr = MPI_Allreduce(&libm_area, ibm_Area,1,MPI_DOUBLE,MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
2724
2725 if (NumberOfBodies > 1) {
2726 ierr = MPI_Allreduce(lIB_Flux,IB_Flux,NumberOfBodies,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD); CHKERRMPI(ierr);
2727 ierr = MPI_Allreduce(lIB_area,IB_Area,NumberOfBodies,MPI_DOUBLE,MPI_SUM,PETSC_COMM_WORLD); CHKERRMPI(ierr);
2728 }
2729
2730 PetscReal correction;
2731
2732 PetscReal *Correction = NULL;
2733 if (NumberOfBodies > 1) {
2734 Correction=(PetscReal *)calloc(NumberOfBodies,sizeof(PetscReal));
2735 for (ibi=0; ibi<NumberOfBodies; ibi++) Correction[ibi]=0.0;
2736 }
2737
2738 if (*ibm_Area > 1.e-15) {
2739 if (flg>1)
2740 correction = (*ibm_Flux + user->
FluxIntpSum)/ ibm_Flux_abs;
2741 else if (flg)
2742 correction = (*ibm_Flux + user->
FluxIntpSum) / *ibm_Area;
2743 else
2744 correction = *ibm_Flux / *ibm_Area;
2745 if (NumberOfBodies > 1)
2746 for (ibi=0; ibi<NumberOfBodies; ibi++) if (IB_Area[ibi]>1.e-15) Correction[ibi] = IB_Flux[ibi] / IB_Area[ibi];
2747 }
2748 else {
2749 correction = 0;
2750 }
2751
2752 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"IBM Uncorrected Flux: %g, Area: %g, Correction: %g\n", *ibm_Flux, *ibm_Area, correction);
2753 if (NumberOfBodies>1){
2754 for (ibi=0; ibi<NumberOfBodies; ibi++)
LOG_ALLOW(
GLOBAL,
LOG_INFO,
" [Body %d] Uncorrected Flux: %g, Area: %g, Correction: %g\n", ibi, IB_Flux[ibi], IB_Area[ibi], Correction[ibi]);
2755 }
2756
2757
2758
2759
2760
2762
2763 for (k=lzs; k<lze; k++) {
2764 for (j=lys; j<lye; j++) {
2765 for (i=lxs; i<lxe; i++) {
2766 if (nvert[k][j][i] < 0.1) {
2767 if (nvert[k][j][i+1] > 0.1 && nvert[k][j][i+1] <ibmval && i < xend) {
2768 if (fabs(ucor[k][j][i].x)>epsilon){
2769 if (flg==3)
2770 ucor[k][j][i].x -=correction*fabs(ucor[k][j][i].x)/
2771 sqrt(csi[k][j][i].x * csi[k][j][i].x +
2772 csi[k][j][i].y * csi[k][j][i].y +
2773 csi[k][j][i].z * csi[k][j][i].z);
2774 else if (flg==2)
2775 ucor[k][j][i].x -=correction*fabs(ucor[k][j][i].x);
2776 else if (NumberOfBodies > 1) {
2777 ibi=(int)((nvert[k][j][i+1]-1.0)*1001);
2778 ucor[k][j][i].x -= sqrt(csi[k][j][i].x * csi[k][j][i].x +
2779 csi[k][j][i].y * csi[k][j][i].y +
2780 csi[k][j][i].z * csi[k][j][i].z) *
2781 Correction[ibi];
2782 }
2783 else
2784 ucor[k][j][i].x -= sqrt(csi[k][j][i].x * csi[k][j][i].x +
2785 csi[k][j][i].y * csi[k][j][i].y +
2786 csi[k][j][i].z * csi[k][j][i].z) *
2787 correction;
2788 }
2789 }
2790 if (nvert[k][j+1][i] > 0.1 && nvert[k][j+1][i] < ibmval && j < yend) {
2791 if (fabs(ucor[k][j][i].y)>epsilon) {
2792 if (flg==3)
2793 ucor[k][j][i].y -=correction*fabs(ucor[k][j][i].y)/
2794 sqrt(eta[k][j][i].x * eta[k][j][i].x +
2795 eta[k][j][i].y * eta[k][j][i].y +
2796 eta[k][j][i].z * eta[k][j][i].z);
2797 else if (flg==2)
2798 ucor[k][j][i].y -=correction*fabs(ucor[k][j][i].y);
2799 else if (NumberOfBodies > 1) {
2800 ibi=(int)((nvert[k][j+1][i]-1.0)*1001);
2801 ucor[k][j][i].y -= sqrt(eta[k][j][i].x * eta[k][j][i].x +
2802 eta[k][j][i].y * eta[k][j][i].y +
2803 eta[k][j][i].z * eta[k][j][i].z) *
2804 Correction[ibi];
2805 }
2806 else
2807 ucor[k][j][i].y -= sqrt(eta[k][j][i].x * eta[k][j][i].x +
2808 eta[k][j][i].y * eta[k][j][i].y +
2809 eta[k][j][i].z * eta[k][j][i].z) *
2810 correction;
2811 }
2812 }
2813 if (nvert[k+1][j][i] > 0.1 && nvert[k+1][j][i] < ibmval && k < zend) {
2814 if (fabs(ucor[k][j][i].z)>epsilon) {
2815 if (flg==3)
2816 ucor[k][j][i].z -= correction*fabs(ucor[k][j][i].z)/
2817 sqrt(zet[k][j][i].x * zet[k][j][i].x +
2818 zet[k][j][i].y * zet[k][j][i].y +
2819 zet[k][j][i].z * zet[k][j][i].z);
2820 else if (flg==2)
2821 ucor[k][j][i].z -= correction*fabs(ucor[k][j][i].z);
2822 else if (NumberOfBodies > 1) {
2823 ibi=(int)((nvert[k+1][j][i]-1.0)*1001);
2824 ucor[k][j][i].z -= sqrt(zet[k][j][i].x * zet[k][j][i].x +
2825 zet[k][j][i].y * zet[k][j][i].y +
2826 zet[k][j][i].z * zet[k][j][i].z) *
2827 Correction[ibi];
2828 }
2829 else
2830 ucor[k][j][i].z -= sqrt(zet[k][j][i].x * zet[k][j][i].x +
2831 zet[k][j][i].y * zet[k][j][i].y +
2832 zet[k][j][i].z * zet[k][j][i].z) *
2833 correction;
2834 }
2835 }
2836 }
2837
2838 if (nvert[k][j][i] > 0.1 && nvert[k][j][i] < ibmval) {
2839 if (nvert[k][j][i+1] < 0.1 && i < xend) {
2840 if (fabs(ucor[k][j][i].x)>epsilon) {
2841 if (flg==3)
2842 ucor[k][j][i].x += correction*fabs(ucor[k][j][i].x)/
2843 sqrt(csi[k][j][i].x * csi[k][j][i].x +
2844 csi[k][j][i].y * csi[k][j][i].y +
2845 csi[k][j][i].z * csi[k][j][i].z);
2846 else if (flg==2)
2847 ucor[k][j][i].x += correction*fabs(ucor[k][j][i].x);
2848 else if (NumberOfBodies > 1) {
2849 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2850 ucor[k][j][i].x += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2851 csi[k][j][i].y * csi[k][j][i].y +
2852 csi[k][j][i].z * csi[k][j][i].z) *
2853 Correction[ibi];
2854 }
2855 else
2856 ucor[k][j][i].x += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2857 csi[k][j][i].y * csi[k][j][i].y +
2858 csi[k][j][i].z * csi[k][j][i].z) *
2859 correction;
2860 }
2861 }
2862 if (nvert[k][j+1][i] < 0.1 && j < yend) {
2863 if (fabs(ucor[k][j][i].y)>epsilon) {
2864 if (flg==3)
2865 ucor[k][j][i].y +=correction*fabs(ucor[k][j][i].y)/
2866 sqrt(eta[k][j][i].x * eta[k][j][i].x +
2867 eta[k][j][i].y * eta[k][j][i].y +
2868 eta[k][j][i].z * eta[k][j][i].z);
2869 else if (flg==2)
2870 ucor[k][j][i].y +=correction*fabs(ucor[k][j][i].y);
2871 else if (NumberOfBodies > 1) {
2872 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2873 ucor[k][j][i].y += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2874 eta[k][j][i].y * eta[k][j][i].y +
2875 eta[k][j][i].z * eta[k][j][i].z) *
2876 Correction[ibi];
2877 }
2878 else
2879 ucor[k][j][i].y += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2880 eta[k][j][i].y * eta[k][j][i].y +
2881 eta[k][j][i].z * eta[k][j][i].z) *
2882 correction;
2883 }
2884 }
2885 if (nvert[k+1][j][i] < 0.1 && k < zend) {
2886 if (fabs(ucor[k][j][i].z)>epsilon) {
2887 if (flg==3)
2888 ucor[k][j][i].z += correction*fabs(ucor[k][j][i].z)/
2889 sqrt(zet[k][j][i].x * zet[k][j][i].x +
2890 zet[k][j][i].y * zet[k][j][i].y +
2891 zet[k][j][i].z * zet[k][j][i].z);
2892 else if (flg==2)
2893 ucor[k][j][i].z += correction*fabs(ucor[k][j][i].z);
2894 else if (NumberOfBodies > 1) {
2895 ibi=(int)((nvert[k][j][i]-1.0)*1001);
2896 ucor[k][j][i].z += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2897 zet[k][j][i].y * zet[k][j][i].y +
2898 zet[k][j][i].z * zet[k][j][i].z) *
2899 Correction[ibi];
2900 }
2901 else
2902 ucor[k][j][i].z += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2903 zet[k][j][i].y * zet[k][j][i].y +
2904 zet[k][j][i].z * zet[k][j][i].z) *
2905 correction;
2906 }
2907 }
2908 }
2909
2910 }
2911 }
2912 }
2913
2914
2915
2916
2917
2919
2920 libm_Flux = 0;
2921 libm_area = 0;
2922 for (k=lzs; k<lze; k++) {
2923 for (j=lys; j<lye; j++) {
2924 for (i=lxs; i<lxe; i++) {
2925 if (nvert[k][j][i] < 0.1) {
2926 if (nvert[k][j][i+1] > 0.1 && nvert[k][j][i+1] < ibmval && i < xend) {
2927 libm_Flux += ucor[k][j][i].x;
2928 libm_area += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2929 csi[k][j][i].y * csi[k][j][i].y +
2930 csi[k][j][i].z * csi[k][j][i].z);
2931
2932 }
2933 if (nvert[k][j+1][i] > 0.1 && nvert[k][j+1][i] < ibmval && j < yend) {
2934 libm_Flux += ucor[k][j][i].y;
2935 libm_area += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2936 eta[k][j][i].y * eta[k][j][i].y +
2937 eta[k][j][i].z * eta[k][j][i].z);
2938 }
2939 if (nvert[k+1][j][i] > 0.1 && nvert[k+1][j][i] < ibmval && k < zend) {
2940 libm_Flux += ucor[k][j][i].z;
2941 libm_area += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2942 zet[k][j][i].y * zet[k][j][i].y +
2943 zet[k][j][i].z * zet[k][j][i].z);
2944 }
2945 }
2946
2947 if (nvert[k][j][i] > 0.1 && nvert[k][j][i] < ibmval) {
2948 if (nvert[k][j][i+1] < 0.1 && i < xend) {
2949 libm_Flux -= ucor[k][j][i].x;
2950 libm_area += sqrt(csi[k][j][i].x * csi[k][j][i].x +
2951 csi[k][j][i].y * csi[k][j][i].y +
2952 csi[k][j][i].z * csi[k][j][i].z);
2953
2954 }
2955 if (nvert[k][j+1][i] < 0.1 && j < yend) {
2956 libm_Flux -= ucor[k][j][i].y;
2957 libm_area += sqrt(eta[k][j][i].x * eta[k][j][i].x +
2958 eta[k][j][i].y * eta[k][j][i].y +
2959 eta[k][j][i].z * eta[k][j][i].z);
2960 }
2961 if (nvert[k+1][j][i] < 0.1 && k < zend) {
2962 libm_Flux -= ucor[k][j][i].z;
2963 libm_area += sqrt(zet[k][j][i].x * zet[k][j][i].x +
2964 zet[k][j][i].y * zet[k][j][i].y +
2965 zet[k][j][i].z * zet[k][j][i].z);
2966 }
2967 }
2968
2969 }
2970 }
2971 }
2972
2973 ierr = MPI_Allreduce(&libm_Flux, ibm_Flux,1,MPI_DOUBLE,MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
2974 ierr = MPI_Allreduce(&libm_area, ibm_Area,1,MPI_DOUBLE,MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
2975
2976
2977
2979
2980
2982 if (xe==mx){
2983 i=mx-2;
2984 for (k=lzs; k<lze; k++) {
2985 for (j=lys; j<lye; j++) {
2986
2987 if ((nvert[k][j][i]>ibmval && nvert[k][j][i+1]<0.1) || (nvert[k][j][i]<0.1 && nvert[k][j][i+1]>ibmval)) ucor[k][j][i].x=0.0;
2988
2989
2990 }
2991 }
2992 }
2993 }
2994
2996 if (ye==my){
2997 j=my-2;
2998 for (k=lzs; k<lze; k++) {
2999 for (i=lxs; i<lxe; i++) {
3000
3001 if ((nvert[k][j][i]>ibmval && nvert[k][j+1][i]<0.1) || (nvert[k][j][i]<0.1 && nvert[k][j+1][i]>ibmval)) ucor[k][j][i].y=0.0;
3002
3003 }
3004 }
3005 }
3006 }
3007
3009 if (ze==mz){
3010 k=mz-2;
3011 for (j=lys; j<lye; j++) {
3012 for (i=lxs; i<lxe; i++) {
3013
3014 if ((nvert[k][j][i]>ibmval && nvert[k+1][j][i]<0.1) || (nvert[k][j][i]<0.1 && nvert[k+1][j][i]>ibmval)) ucor[k][j][i].z=0.0;
3015
3016 }
3017 }
3018 }
3019 }
3020
3021
3022 DMDAVecRestoreArray(da, user->
lNvert, &nvert);
3023 DMDAVecRestoreArray(fda, user->
lCsi, &csi);
3024 DMDAVecRestoreArray(fda, user->
lEta, &eta);
3025 DMDAVecRestoreArray(fda, user->
lZet, &zet);
3026 DMDAVecRestoreArray(fda, user->
Ucont, &ucor);
3027
3028 DMGlobalToLocalBegin(fda, user->
Ucont, INSERT_VALUES, user->
lUcont);
3029 DMGlobalToLocalEnd(fda, user->
Ucont, INSERT_VALUES, user->
lUcont);
3030
3031
3032
3033 PetscBool has_periodic_ucont_seam =
3040 PetscInt periodic_seam_updates_local = 0, periodic_seam_updates_global = 0;
3041 PetscReal periodic_seam_max_delta_local = 0.0, periodic_seam_max_delta_global = 0.0;
3042
3043 if (has_periodic_ucont_seam) {
3045 }
3046
3047 DMDAVecGetArray(fda, user->
lUcont, &lucor);
3048 DMDAVecGetArray(fda, user->
Ucont, &ucor);
3049
3051 if (xs==0){
3052 i=xs;
3053 for (k=zs; k<ze; k++) {
3054 for (j=ys; j<ye; j++) {
3055 if(j>0 && k>0 && j<user->JM && k<user->KM){
3056 PetscReal old_val = ucor[k][j][i].x;
3057 PetscReal new_val = lucor[k][j][i-2].x;
3058 PetscReal delta = PetscAbsReal(new_val - old_val);
3059 ucor[k][j][i].x = new_val;
3060 if (delta > 1.0e-14) periodic_seam_updates_local++;
3061 periodic_seam_max_delta_local = PetscMax(periodic_seam_max_delta_local, delta);
3062 }
3063 }
3064 }
3065 }
3066 }
3068 if (ys==0){
3069 j=ys;
3070 for (k=zs; k<ze; k++) {
3071 for (i=xs; i<xe; i++) {
3072 if(i>0 && k>0 && i<user->IM && k<user->KM){
3073 PetscReal old_val = ucor[k][j][i].y;
3074 PetscReal new_val = lucor[k][j-2][i].y;
3075 PetscReal delta = PetscAbsReal(new_val - old_val);
3076 ucor[k][j][i].y = new_val;
3077 if (delta > 1.0e-14) periodic_seam_updates_local++;
3078 periodic_seam_max_delta_local = PetscMax(periodic_seam_max_delta_local, delta);
3079 }
3080 }
3081 }
3082 }
3083 }
3085 if (zs==0){
3086 k=zs;
3087 for (j=ys; j<ye; j++) {
3088 for (i=xs; i<xe; i++) {
3089 if(i>0 && j>0 && i<user->IM && j<user->JM){
3090 PetscReal old_val = ucor[k][j][i].z;
3091 PetscReal new_val = lucor[k-2][j][i].z;
3092 PetscReal delta = PetscAbsReal(new_val - old_val);
3093 ucor[k][j][i].z = new_val;
3094 if (delta > 1.0e-14) periodic_seam_updates_local++;
3095 periodic_seam_max_delta_local = PetscMax(periodic_seam_max_delta_local, delta);
3096 }
3097 }
3098 }
3099 }
3100 }
3101 DMDAVecRestoreArray(fda, user->
lUcont, &lucor);
3102 DMDAVecRestoreArray(fda, user->
Ucont, &ucor);
3103
3104 if (has_periodic_ucont_seam) {
3105 ierr = MPI_Allreduce(&periodic_seam_updates_local, &periodic_seam_updates_global, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD); CHKERRMPI(ierr);
3106 ierr = MPI_Allreduce(&periodic_seam_max_delta_local, &periodic_seam_max_delta_global, 1, MPIU_REAL, MPI_MAX, PETSC_COMM_WORLD); CHKERRMPI(ierr);
3108 "VolumeFlux: legacy periodic Ucont seam refresh updated %d values, max |delta|=%.6e.\n",
3109 (int)periodic_seam_updates_global, periodic_seam_max_delta_global);
3110 }
3111
3112
3113
3114 DMGlobalToLocalBegin(user->
fda, user->
Ucont, INSERT_VALUES, user->
lUcont);
3115 DMGlobalToLocalEnd(user->
fda, user->
Ucont, INSERT_VALUES, user->
lUcont);
3116
3117 if (NumberOfBodies > 1) {
3118 free(lIB_Flux);
3119 free(lIB_area);
3120 free(IB_Flux);
3121 free(IB_Area);
3122 free(Correction);
3123 }
3124
3126
3127 return 0;
3128}