Skip to content

Commit

Permalink
Fortran Interfaces: Add new average down functions (#4124)
Browse files Browse the repository at this point in the history
Add average down function for cell-centered data without volume
weighting.

Add average down function for nodal data.
  • Loading branch information
WeiqunZhang committed Sep 4, 2024
1 parent dea9bb1 commit 65d10a1
Show file tree
Hide file tree
Showing 2 changed files with 33 additions and 3 deletions.
10 changes: 8 additions & 2 deletions Src/F_Interfaces/Base/AMReX_multifabutil_fi.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -5,12 +5,18 @@ using namespace amrex;
extern "C"
{
void amrex_fi_average_down (const MultiFab* S_fine, MultiFab* S_crse,
const Geometry* fgeom, const Geometry* cgeom,
int scomp, int ncomp, int rr)
const Geometry* fgeom, const Geometry* cgeom,
int scomp, int ncomp, int rr)
{
amrex::average_down(*S_fine, *S_crse, *fgeom, *cgeom, scomp, ncomp, rr);
}

void amrex_fi_average_down_cell_node (const MultiFab* S_fine, MultiFab* S_crse,
int scomp, int ncomp, int rr)
{
amrex::average_down(*S_fine, *S_crse, scomp, ncomp, rr);
}

void amrex_fi_average_down_faces (MultiFab const* fmf[], MultiFab* cmf[],
Geometry const* cgeom, int scomp, int ncomp,
int rr)
Expand Down
26 changes: 25 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,11 @@ module amrex_multifabutil_module
implicit none
private

public :: amrex_average_down, amrex_average_down_faces, amrex_average_cellcenter_to_face
public :: amrex_average_down, & ! volume weighted average down of cell data
& amrex_average_down_cell, & ! average down of cell data
& amrex_average_down_node, & ! average down of nodal data
& amrex_average_down_faces, & ! average down of face data
& amrex_average_cellcenter_to_face ! average from cell centers to faces

interface
subroutine amrex_fi_average_down (fmf, cmf, fgeom, cgeom, scomp, ncomp, rr) bind(c)
Expand All @@ -18,6 +22,13 @@ 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_cell_node (fmf, cmf, scomp, ncomp, rr) bind(c)
import
implicit none
type(c_ptr), value :: fmf, cmf
integer(c_int), value :: scomp, ncomp, rr
end subroutine amrex_fi_average_down_cell_node

subroutine amrex_fi_average_down_faces (fmf, cmf, cgeom, scomp, ncomp, rr) bind(c)
import
implicit none
Expand Down Expand Up @@ -45,6 +56,19 @@ subroutine amrex_average_down (fmf, cmf, fgeom, cgeom, scomp, ncomp, rr)
call amrex_fi_average_down(fmf%p, cmf%p, fgeom%p, cgeom%p, scomp-1, ncomp, rr)
end subroutine amrex_average_down

subroutine amrex_average_down_cell (fmf, cmf, scomp, ncomp, rr)
type(amrex_multifab), intent(in ) :: fmf
type(amrex_multifab), intent(inout) :: cmf
integer, intent(in) :: scomp, ncomp, rr
call amrex_fi_average_down_cell_node(fmf%p, cmf%p, scomp-1, ncomp, rr)
end subroutine amrex_average_down_cell

subroutine amrex_average_down_node (fmf, cmf, scomp, ncomp, rr)
type(amrex_multifab), intent(in ) :: fmf
type(amrex_multifab), intent(inout) :: cmf
integer, intent(in) :: scomp, ncomp, rr
call amrex_fi_average_down_cell_node(fmf%p, cmf%p, scomp-1, ncomp, rr)
end subroutine amrex_average_down_node

subroutine amrex_average_down_faces (fmf, cmf, cgeom, scomp, ncomp, rr)
type(amrex_multifab), intent(in ) :: fmf(amrex_spacedim)
Expand Down

0 comments on commit 65d10a1

Please sign in to comment.