Skip to content

Commit

Permalink
Add Fortran interface for average_down_faces (AMReX-Codes#3553)
Browse files Browse the repository at this point in the history
  • Loading branch information
WeiqunZhang authored Sep 25, 2023
1 parent 2e99628 commit f08b40b
Show file tree
Hide file tree
Showing 3 changed files with 42 additions and 2 deletions.
3 changes: 2 additions & 1 deletion Src/Base/AMReX_MultiFabUtil.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -335,7 +335,8 @@ namespace amrex
amrex::average_down(S_fine, S_crse, scomp, ncomp, ratio);
#else

AMREX_ASSERT(S_crse.nComp() == S_fine.nComp());
AMREX_ASSERT(S_crse.nComp() == S_fine.nComp() &&
S_fine.is_cell_centered() && S_crse.is_cell_centered());

//
// Coarsen() the fine stuff on processors owning the fine data.
Expand Down
16 changes: 16 additions & 0 deletions Src/F_Interfaces/Base/AMReX_multifabutil_fi.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,22 @@ extern "C"
amrex::average_down(*S_fine, *S_crse, *fgeom, *cgeom, scomp, ncomp, rr);
}

void amrex_fi_average_down_faces (MultiFab const* fmf[], MultiFab* cmf[],
Geometry const* cgeom, int scomp, int ncomp,
int rr)
{
Array<MultiFab,AMREX_SPACEDIM> fine
{AMREX_D_DECL(MultiFab(*fmf[0], amrex::make_alias, scomp, ncomp),
MultiFab(*fmf[1], amrex::make_alias, scomp, ncomp),
MultiFab(*fmf[2], amrex::make_alias, scomp, ncomp))};
Array<MultiFab,AMREX_SPACEDIM> crse
{AMREX_D_DECL(MultiFab(*cmf[0], amrex::make_alias, scomp, ncomp),
MultiFab(*cmf[1], amrex::make_alias, scomp, ncomp),
MultiFab(*cmf[2], amrex::make_alias, scomp, ncomp))};
amrex::average_down_faces(GetArrOfConstPtrs(fine), GetArrOfPtrs(crse),
IntVect(rr), *cgeom);
}

void amrex_fi_average_cellcenter_to_face (MultiFab* fc[], const MultiFab* cc, const Geometry* geom)
{
Vector<MultiFab*> fcv {AMREX_D_DECL(fc[0], fc[1], fc[2])};
Expand Down
25 changes: 24 additions & 1 deletion Src/F_Interfaces/Base/AMReX_multifabutil_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@ module amrex_multifabutil_module
implicit none
private

public :: amrex_average_down, amrex_average_cellcenter_to_face
public :: amrex_average_down, amrex_average_down_faces, amrex_average_cellcenter_to_face

interface
subroutine amrex_fi_average_down (fmf, cmf, fgeom, cgeom, scomp, ncomp, rr) bind(c)
Expand All @@ -18,6 +18,14 @@ subroutine amrex_fi_average_down (fmf, cmf, fgeom, cgeom, scomp, ncomp, rr) bind
integer(c_int), value :: scomp, ncomp, rr
end subroutine amrex_fi_average_down

subroutine amrex_fi_average_down_faces (fmf, cmf, cgeom, scomp, ncomp, rr) bind(c)
import
implicit none
type(c_ptr), intent(in) :: fmf(*), cmf(*)
type(c_ptr), value :: cgeom
integer(c_int), value :: scomp, ncomp, rr
end subroutine amrex_fi_average_down_faces

subroutine amrex_fi_average_cellcenter_to_face (facemf, ccmf, geom) bind(c)
import
implicit none
Expand All @@ -38,6 +46,21 @@ subroutine amrex_average_down (fmf, cmf, fgeom, cgeom, scomp, ncomp, rr)
end subroutine amrex_average_down


subroutine amrex_average_down_faces (fmf, cmf, cgeom, scomp, ncomp, rr)
type(amrex_multifab), intent(in ) :: fmf(amrex_spacedim)
type(amrex_multifab), intent(inout) :: cmf(amrex_spacedim)
type(amrex_geometry), intent(in) :: cgeom ! coarse level geometry
integer, intent(in) :: scomp, ncomp, rr
type(c_ptr) :: cp(amrex_spacedim), fp(amrex_spacedim)
integer :: idim
do idim = 1, amrex_spacedim
fp(idim) = fmf(idim)%p
cp(idim) = cmf(idim)%p
end do
call amrex_fi_average_down_faces(fp, cp, cgeom%p, scomp-1, ncomp, rr)
end subroutine amrex_average_down_faces


subroutine amrex_average_cellcenter_to_face (facemf, ccmf, geom)
type(amrex_multifab), intent(inout) :: facemf(amrex_spacedim)
type(amrex_multifab), intent(in) :: ccmf
Expand Down

0 comments on commit f08b40b

Please sign in to comment.