From 83497899f0a2b4b8d3b1ff5c1ded0d22ad85a749 Mon Sep 17 00:00:00 2001 From: Sam Dunham <9784969+dunhamsj@users.noreply.github.com> Date: Tue, 1 Oct 2024 15:27:49 -0400 Subject: [PATCH] Add comparison operator for boxarray and distromap. Add hdf5 to dep.py (#4173) Add operator(==) for boxarray and for distromap for the Fortran interface. Also add "hdf5" to dep.py to suppress warnings. --------- Co-authored-by: Weiqun Zhang --- Src/F_Interfaces/Base/AMReX_boxarray_fi.cpp | 5 +++++ Src/F_Interfaces/Base/AMReX_boxarray_mod.F90 | 21 +++++++++++++++++-- Src/F_Interfaces/Base/AMReX_distromap_fi.cpp | 5 +++++ Src/F_Interfaces/Base/AMReX_distromap_mod.F90 | 18 +++++++++++++++- Tools/F_scripts/dep.py | 2 +- 5 files changed, 47 insertions(+), 4 deletions(-) diff --git a/Src/F_Interfaces/Base/AMReX_boxarray_fi.cpp b/Src/F_Interfaces/Base/AMReX_boxarray_fi.cpp index dd7916a9adf..48d8eae7a6f 100644 --- a/Src/F_Interfaces/Base/AMReX_boxarray_fi.cpp +++ b/Src/F_Interfaces/Base/AMReX_boxarray_fi.cpp @@ -86,4 +86,9 @@ extern "C" { Box bx(IntVect(lo), IntVect(hi), ba->ixType()); return ba->intersects(bx); } + + int amrex_fi_boxarray_issame (const BoxArray* baa, const BoxArray* bab) + { + return *baa == *bab; + } } diff --git a/Src/F_Interfaces/Base/AMReX_boxarray_mod.F90 b/Src/F_Interfaces/Base/AMReX_boxarray_mod.F90 index 0181c6cfb9c..f85b5e8d74e 100644 --- a/Src/F_Interfaces/Base/AMReX_boxarray_mod.F90 +++ b/Src/F_Interfaces/Base/AMReX_boxarray_mod.F90 @@ -9,7 +9,8 @@ module amrex_boxarray_module private - public :: amrex_boxarray_build, amrex_boxarray_destroy, amrex_print + public :: amrex_boxarray_build, amrex_boxarray_destroy, amrex_print, & + operator(==) type, public :: amrex_boxarray logical :: owner = .false. @@ -36,6 +37,10 @@ module amrex_boxarray_module #endif end type amrex_boxarray + interface operator(==) + module procedure amrex_boxarray_issame + end interface operator(==) + interface amrex_boxarray_build module procedure amrex_boxarray_build_bx module procedure amrex_boxarray_build_bxs @@ -122,12 +127,18 @@ pure function amrex_fi_boxarray_numpts (ba) bind(c) integer(amrex_long) :: amrex_fi_boxarray_numpts end function amrex_fi_boxarray_numpts - pure integer function amrex_fi_boxarray_intersects_box (ba, lo, hi) bind(c) + pure integer(c_int) function amrex_fi_boxarray_intersects_box (ba, lo, hi) bind(c) import implicit none type(c_ptr), value, intent(in) :: ba integer, intent(in) :: lo(*), hi(*) end function amrex_fi_boxarray_intersects_box + + pure integer(c_int) function amrex_fi_boxarray_issame (baa, bab) bind(c) + import + implicit none + type(c_ptr), value, intent(in) :: baa, bab + end function amrex_fi_boxarray_issame end interface contains @@ -258,4 +269,10 @@ pure function amrex_boxarray_intersects_box (this, bx) result(r) r = ir .ne. 0 end function amrex_boxarray_intersects_box + pure logical function amrex_boxarray_issame(baa, bab) result(r) + type(amrex_boxarray), intent(in) :: baa + type(amrex_boxarray), intent(in) :: bab + r = amrex_fi_boxarray_issame(baa%p, bab%p) .ne. 0 + end function amrex_boxarray_issame + end module amrex_boxarray_module diff --git a/Src/F_Interfaces/Base/AMReX_distromap_fi.cpp b/Src/F_Interfaces/Base/AMReX_distromap_fi.cpp index e50031a5887..7fc7adc1718 100644 --- a/Src/F_Interfaces/Base/AMReX_distromap_fi.cpp +++ b/Src/F_Interfaces/Base/AMReX_distromap_fi.cpp @@ -41,4 +41,9 @@ extern "C" { { AllPrint() << *dm; } + + int amrex_fi_distromap_issame (const DistributionMapping* dma, const DistributionMapping* dmb) + { + return *dma == *dmb; + } } diff --git a/Src/F_Interfaces/Base/AMReX_distromap_mod.F90 b/Src/F_Interfaces/Base/AMReX_distromap_mod.F90 index adbb91b4421..9c0884168ea 100644 --- a/Src/F_Interfaces/Base/AMReX_distromap_mod.F90 +++ b/Src/F_Interfaces/Base/AMReX_distromap_mod.F90 @@ -8,7 +8,8 @@ module amrex_distromap_module private - public :: amrex_distromap_build, amrex_distromap_destroy, amrex_print + public :: amrex_distromap_build, amrex_distromap_destroy, amrex_print, & + operator(==) type, public :: amrex_distromap logical :: owner = .false. @@ -25,6 +26,10 @@ module amrex_distromap_module #endif end type amrex_distromap + interface operator(==) + module procedure amrex_distromap_issame + end interface operator(==) + interface amrex_distromap_build module procedure amrex_distromap_build_ba module procedure amrex_distromap_build_pmap @@ -89,6 +94,12 @@ subroutine amrex_fi_print_distromap (dm) bind(c) implicit none type(c_ptr), value :: dm end subroutine amrex_fi_print_distromap + + pure integer(c_int) function amrex_fi_distromap_issame (dma, dmb) bind(c) + import + implicit none + type(c_ptr), value, intent(in) :: dma, dmb + end function amrex_fi_distromap_issame end interface contains @@ -158,4 +169,9 @@ subroutine amrex_distromap_print (dm) call amrex_fi_print_distromap(dm%p) end subroutine amrex_distromap_print + pure logical function amrex_distromap_issame (dma, dmb) result(r) + type(amrex_distromap), intent(in) :: dma, dmb + r = amrex_fi_distromap_issame(dma%p, dmb%p) .ne. 0 + end function amrex_distromap_issame + end module amrex_distromap_module diff --git a/Tools/F_scripts/dep.py b/Tools/F_scripts/dep.py index 151d1a9f9ab..e5eb74cb40b 100755 --- a/Tools/F_scripts/dep.py +++ b/Tools/F_scripts/dep.py @@ -28,7 +28,7 @@ import preprocess # modules to ignore in the dependencies -IGNORES = ["iso_c_binding", "iso_fortran_env", "omp_lib", "mpi", "cudafor", "openacc", "hdf"] +IGNORES = ["iso_c_binding", "iso_fortran_env", "omp_lib", "mpi", "cudafor", "openacc", "hdf", "hdf5"] # regular expression for "{}module{}name", where {} can be any number # of spaces. We use 4 groups here, denoted by (), so the name of the