diff --git a/amr-wind/utilities/sampling/ProbeSampler.H b/amr-wind/utilities/sampling/ProbeSampler.H index de59322241..c432563ca9 100644 --- a/amr-wind/utilities/sampling/ProbeSampler.H +++ b/amr-wind/utilities/sampling/ProbeSampler.H @@ -52,6 +52,9 @@ private: const CFDSim& m_sim; SampleLocType m_probes; + amrex::Vector m_offset_vector{0.0, 0.0, 0.0}; + amrex::Vector m_poffsets; + std::string m_label; int m_id{-1}; int m_npts{0}; diff --git a/amr-wind/utilities/sampling/ProbeSampler.cpp b/amr-wind/utilities/sampling/ProbeSampler.cpp index dfe7aed42e..7d033f6268 100644 --- a/amr-wind/utilities/sampling/ProbeSampler.cpp +++ b/amr-wind/utilities/sampling/ProbeSampler.cpp @@ -20,13 +20,38 @@ void ProbeSampler::initialize(const std::string& key) amrex::Abort("Cannot find probe location file: " + pfile); } - ifh >> m_npts; + pp.queryarr("offsets", m_poffsets); + if (m_poffsets.size() > 0) { + pp.getarr("offset_vector", m_offset_vector); + AMREX_ALWAYS_ASSERT( + static_cast(m_offset_vector.size()) == AMREX_SPACEDIM); + } else { + // No offsets is implemented as 1 offset of 0. + m_poffsets.push_back(0.0); + } + + int npts_file = 0; + ifh >> npts_file; ifh.ignore(std::numeric_limits::max(), '\n'); + SampleLocType probes_file; + probes_file.resize(npts_file); + m_npts = m_poffsets.size() * npts_file; m_probes.resize(m_npts); - for (int i = 0; i < m_npts; ++i) { - ifh >> m_probes[i][0] >> m_probes[i][1] >> m_probes[i][2]; + // Read through points in file + for (int i = 0; i < npts_file; ++i) { + ifh >> probes_file[i][0] >> probes_file[i][1] >> probes_file[i][2]; ifh.ignore(std::numeric_limits::max(), '\n'); } + // Incorporate offsets + for (int n = 0; n < m_poffsets.size(); ++n) { + for (int i = 0; i < npts_file; ++i) { + for (int d = 0; d < AMREX_SPACEDIM; ++d) { + m_probes[i + n * npts_file][d] = + probes_file[i][d] + m_poffsets[n] * m_offset_vector[d]; + } + } + } + check_bounds(); } @@ -69,6 +94,8 @@ void ProbeSampler::sampling_locations(SampleLocType& locs) const void ProbeSampler::define_netcdf_metadata(const ncutils::NCGroup& grp) const { grp.put_attr("sampling_type", identifier()); + grp.put_attr("offset_vector", m_offset_vector); + grp.put_attr("offsets", m_poffsets); } #else void ProbeSampler::define_netcdf_metadata( diff --git a/docs/sphinx/user/inputs_Sampling.rst b/docs/sphinx/user/inputs_Sampling.rst index 110d2252b4..e26eec32fc 100644 --- a/docs/sphinx/user/inputs_Sampling.rst +++ b/docs/sphinx/user/inputs_Sampling.rst @@ -158,6 +158,10 @@ Example:: The first line of the file contains the total number of probes for this set. This is followed by the coordinates (three real numbers), one line per probe. +This type of sampler also supports the ``offset_vector`` and ``offsets`` options +implemented with the plane sampler, shown above. For the probe sampler, +these options apply offsets to the positions of all the points provided in the +probe location file. Sampling on a volume ````````````````````` diff --git a/unit_tests/utilities/test_sampling.cpp b/unit_tests/utilities/test_sampling.cpp index dbc453a910..5f12a73370 100644 --- a/unit_tests/utilities/test_sampling.cpp +++ b/unit_tests/utilities/test_sampling.cpp @@ -3,6 +3,7 @@ #include "amr-wind/utilities/sampling/Sampling.H" #include "amr-wind/utilities/sampling/SamplingContainer.H" +#include "amr-wind/utilities/sampling/ProbeSampler.H" #include "amr-wind/utilities/sampling/PlaneSampler.H" #include "amr-wind/utilities/sampling/VolumeSampler.H" #include "amr-wind/utilities/sampling/DTUSpinnerSampler.H" @@ -61,6 +62,17 @@ void init_int_field(amr_wind::IntField& fld) amrex::Gpu::synchronize(); } +void write_probe_sampler_file(const std::string& fname) +{ + std::ofstream os(fname); + // Total number of points + os << "3\n"; + // Coordinates + os << "0.0\t0.0\t0.0\n"; + os << "60.0\t2.0\t3.0\n"; + os << "100.0\t8.0\t5.0\n"; +} + class SamplingImpl : public amr_wind::sampling::Sampling { public: @@ -264,6 +276,52 @@ TEST_F(SamplingTest, sampling_timing) EXPECT_TRUE(probes.write_flag); } +TEST_F(SamplingTest, probe_sampler) +{ + initialize_mesh(); + + constexpr amrex::Real tol = 1.0e-12; + std::string fname = "probes.txt"; + // Write file + write_probe_sampler_file(fname); + + { + amrex::ParmParse pp("cloud"); + pp.add("probe_location_file", fname); + pp.addarr("offsets", amrex::Vector{1.0, 2.5}); + pp.addarr("offset_vector", amrex::Vector{0.2, 0.5, 1.0}); + } + + amr_wind::sampling::ProbeSampler cloud(sim()); + cloud.initialize("cloud"); + amr_wind::sampling::ProbeSampler::SampleLocType locs; + cloud.sampling_locations(locs); + + ASSERT_EQ(locs.size(), 3 * 2); + const amrex::Vector xprobe_golds{0.2, 60.2, 100.2, + 0.5, 60.5, 100.5}; + const amrex::Vector yprobe_golds{0.5, 2.5, 8.5, + 1.25, 3.25, 9.25}; + const amrex::Vector zprobe_golds{1.0, 4.0, 6.0, 2.5, 5.5, 7.5}; + for (int n = 0; n < locs.size(); ++n) { + EXPECT_NEAR(locs[n][0], xprobe_golds[n], tol); + EXPECT_NEAR(locs[n][1], yprobe_golds[n], tol); + EXPECT_NEAR(locs[n][2], zprobe_golds[n], tol); + } + + // Remove file + const char* fname_char = fname.c_str(); + { + std::ifstream f(fname_char); + if (f.good()) { + remove(fname_char); + } + // Check that file is removed + std::ifstream ff(fname_char); + EXPECT_FALSE(ff.good()); + } +} + TEST_F(SamplingTest, plane_sampler) { initialize_mesh();