From: Rene Gassmoeller Date: Tue, 8 Aug 2017 22:02:34 +0000 (-0600) Subject: Address initial comments. Add tests. X-Git-Tag: v9.0.0-rc1~1013^2~5 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=6861db8263a1c69ac795e67215b9c2dd655a16f4;p=dealii.git Address initial comments. Add tests. --- diff --git a/include/deal.II/particle/particle.h b/include/deal.II/particles/particle.h similarity index 88% rename from include/deal.II/particle/particle.h rename to include/deal.II/particles/particle.h index ab488769b3..6c7bc369fd 100644 --- a/include/deal.II/particle/particle.h +++ b/include/deal.II/particles/particle.h @@ -13,10 +13,10 @@ // // --------------------------------------------------------------------- -#ifndef dealii__particle_particle_h -#define dealii__particle_particle_h +#ifndef dealii__particles_particle_h +#define dealii__particles_particle_h -#include +#include #include #include @@ -26,7 +26,11 @@ DEAL_II_NAMESPACE_OPEN -namespace Particle +/** + * A namespace that contains all classes that are related to the particle + * implementation, in particular the fundamental Particle class. + */ +namespace Particles { /** * A namespace for all type definitions related to particles. @@ -53,7 +57,7 @@ namespace Particle * An identifier that denotes the MPI type associated with * types::global_dof_index. */ -# define ASPECT_PARTICLE_INDEX_MPI_TYPE MPI_UNSIGNED_LONG_LONG +# define PARTICLE_INDEX_MPI_TYPE MPI_UNSIGNED_LONG_LONG #else /** * The type used for indices of particles. While in @@ -72,7 +76,7 @@ namespace Particle * An identifier that denotes the MPI type associated with * types::global_dof_index. */ -# define ASPECT_PARTICLE_INDEX_MPI_TYPE MPI_UNSIGNED +# define PARTICLE_INDEX_MPI_TYPE MPI_UNSIGNED #endif } @@ -98,14 +102,14 @@ namespace Particle /** * Constructor for Particle, creates a particle with the specified - * ID at the specified location. Note that Aspect - * does not check for duplicate particle IDs so the generator must + * ID at the specified location. Note that deal.II + * does not check for duplicate particle IDs so the user must * make sure the IDs are unique over all processes. * - * @param[in] new_location Initial location of particle. - * @param[in] new_reference_location Initial location of the particle + * @param new_location Initial location of particle. + * @param new_reference_location Initial location of the particle * in the coordinate system of the reference cell. - * @param[in] new_id Globally unique ID number of particle. + * @param new_id Globally unique ID number of particle. */ Particle (const Point &new_location, const Point &new_reference_location, @@ -195,8 +199,7 @@ namespace Particle get_location () const; /** - * Set the reference location of this particle. Note that this does not - * check whether this is a valid location in the simulation domain. + * Set the reference location of this particle. * * @param [in] new_loc The new reference location for this particle. */ @@ -204,17 +207,13 @@ namespace Particle set_reference_location (const Point &new_loc); /** - * Get the reference location of this particle in its current cell. - * - * @return The reference location of this particle. + * Return the reference location of this particle in its current cell. */ const Point & get_reference_location () const; /** - * Get the ID number of this particle. - * - * @return The id of this particle. + * Return the ID number of this particle. */ types::particle_index get_id () const; @@ -268,8 +267,6 @@ namespace Particle /** * Current particle location in the reference cell. - * Storing this reduces the number of times we need to compute this - * location, which takes a significant amount of computing time. */ Point reference_location; diff --git a/include/deal.II/particle/property_pool.h b/include/deal.II/particles/property_pool.h similarity index 95% rename from include/deal.II/particle/property_pool.h rename to include/deal.II/particles/property_pool.h index 99dc52c028..a392c23449 100644 --- a/include/deal.II/particle/property_pool.h +++ b/include/deal.II/particles/property_pool.h @@ -13,14 +13,14 @@ // // --------------------------------------------------------------------- -#ifndef dealii__particle_property_pool_h -#define dealii__particle_property_pool_h +#ifndef dealii__particles_property_pool_h +#define dealii__particles_property_pool_h #include DEAL_II_NAMESPACE_OPEN -namespace Particle +namespace Particles { /** * This class manages the memory space in which particles store their diff --git a/source/CMakeLists.txt b/source/CMakeLists.txt index e21a353811..d9c364a350 100644 --- a/source/CMakeLists.txt +++ b/source/CMakeLists.txt @@ -48,7 +48,7 @@ ADD_SUBDIRECTORY(integrators) ADD_SUBDIRECTORY(matrix_free) ADD_SUBDIRECTORY(meshworker) ADD_SUBDIRECTORY(opencascade) -ADD_SUBDIRECTORY(particle) +ADD_SUBDIRECTORY(particles) ADD_SUBDIRECTORY(physics) ADD_SUBDIRECTORY(non_matching) ADD_SUBDIRECTORY(sundials) diff --git a/source/particle/CMakeLists.txt b/source/particles/CMakeLists.txt similarity index 90% rename from source/particle/CMakeLists.txt rename to source/particles/CMakeLists.txt index 99db585c88..b45a281df6 100644 --- a/source/particle/CMakeLists.txt +++ b/source/particles/CMakeLists.txt @@ -1,6 +1,6 @@ ## --------------------------------------------------------------------- ## -## Copyright (C) 2016 by the deal.II authors +## Copyright (C) 2017 by the deal.II authors ## ## This file is part of the deal.II library. ## @@ -26,7 +26,7 @@ SET(_inst ) FILE(GLOB _header - ${CMAKE_SOURCE_DIR}/include/deal.II/particle/*.h + ${CMAKE_SOURCE_DIR}/include/deal.II/particles/*.h ) DEAL_II_ADD_LIBRARY(obj_particle OBJECT ${_src} ${_header} ${_inst}) diff --git a/source/particle/particle.cc b/source/particles/particle.cc similarity index 98% rename from source/particle/particle.cc rename to source/particles/particle.cc index f842b41148..11e74a23ca 100644 --- a/source/particle/particle.cc +++ b/source/particles/particle.cc @@ -13,11 +13,11 @@ // // --------------------------------------------------------------------- -#include +#include DEAL_II_NAMESPACE_OPEN -namespace Particle +namespace Particles { template Particle::Particle () @@ -92,8 +92,6 @@ namespace Particle data = static_cast (pdata); } -#ifdef DEAL_II_WITH_CXX11 - template Particle::Particle (Particle &&particle) : @@ -149,7 +147,6 @@ namespace Particle } return *this; } -#endif template Particle::~Particle () diff --git a/source/particle/particle.inst.in b/source/particles/particle.inst.in similarity index 96% rename from source/particle/particle.inst.in rename to source/particles/particle.inst.in index 0c888d702c..4b693348c0 100644 --- a/source/particle/particle.inst.in +++ b/source/particles/particle.inst.in @@ -16,7 +16,7 @@ for (deal_II_dimension : DIMENSIONS) { - namespace Particle + namespace Particles \{ template class Particle ; diff --git a/source/particle/property_pool.cc b/source/particles/property_pool.cc similarity index 88% rename from source/particle/property_pool.cc rename to source/particles/property_pool.cc index 1aa1adaa05..4a85e3fffd 100644 --- a/source/particle/property_pool.cc +++ b/source/particles/property_pool.cc @@ -14,14 +14,13 @@ // --------------------------------------------------------------------- -#include -#include +#include DEAL_II_NAMESPACE_OPEN -namespace Particle +namespace Particles { - const PropertyPool::Handle PropertyPool::invalid_handle = NULL; + const typename PropertyPool::Handle PropertyPool::invalid_handle = NULL; PropertyPool::PropertyPool (const unsigned int n_properties_per_slot) @@ -54,6 +53,7 @@ namespace Particle } + void PropertyPool::reserve(const std::size_t size) { diff --git a/tests/particle/CMakeLists.txt b/tests/particles/CMakeLists.txt similarity index 100% rename from tests/particle/CMakeLists.txt rename to tests/particles/CMakeLists.txt diff --git a/tests/particle/particle_01.cc b/tests/particles/particle_01.cc similarity index 93% rename from tests/particle/particle_01.cc rename to tests/particles/particle_01.cc index fa8a009385..876e03ac16 100644 --- a/tests/particle/particle_01.cc +++ b/tests/particles/particle_01.cc @@ -18,7 +18,7 @@ // check the creation and destruction of particles #include "../tests.h" -#include +#include #include #include @@ -27,7 +27,7 @@ template void test () { { - Particle::Particle particle; + Particles::Particle particle; deallog << "Particle location: " << particle.get_location() << std::endl; diff --git a/tests/particle/particle_01.output b/tests/particles/particle_01.output similarity index 100% rename from tests/particle/particle_01.output rename to tests/particles/particle_01.output diff --git a/tests/particles/particle_02.cc b/tests/particles/particle_02.cc new file mode 100644 index 0000000000..95cc6929f8 --- /dev/null +++ b/tests/particles/particle_02.cc @@ -0,0 +1,72 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// check Particle constructors, copy, and move operations. + +#include "../tests.h" +#include +#include +#include + + +template +void test () +{ + { + Point<2> position; + position(0) = 0.3; + position(1) = 0.5; + + Point<2> reference_position; + reference_position(0) = 0.2; + reference_position(1) = 0.4; + + const Particles::types::particle_index index(7); + + Particles::Particle<2> particle(position,reference_position,index); + + deallog << "Particle location: " << particle.get_location() << std::endl + << "Particle reference location: " << particle.get_reference_location() << std::endl + << "Particle index: " << particle.get_id() << std::endl; + + const Particles::Particle<2> copy(particle); + + deallog << "Copy particle location: " << copy.get_location() << std::endl + << "Copy particle reference location: " << copy.get_reference_location() << std::endl + << "Copy particle index: " << copy.get_id() << std::endl; + + const Particles::Particle<2> moved_particle(std::move(particle)); + + deallog << "Moved particle location: " << moved_particle.get_location() << std::endl + << "Moved particle reference location: " << moved_particle.get_reference_location() << std::endl + << "Moved particle index: " << moved_particle.get_id() << std::endl; + + deallog << "Original particle location: " << particle.get_location() << std::endl + << "Original particle reference location: " << particle.get_reference_location() << std::endl + << "Original particle index: " << particle.get_id() << std::endl; + } + + deallog << "OK" << std::endl; +} + + + +int main () +{ + initlog(); + test<2>(); +} diff --git a/tests/particles/particle_02.output b/tests/particles/particle_02.output new file mode 100644 index 0000000000..2a4fefe16d --- /dev/null +++ b/tests/particles/particle_02.output @@ -0,0 +1,14 @@ + +DEAL::Particle location: 0.300000 0.500000 +DEAL::Particle reference location: 0.200000 0.400000 +DEAL::Particle index: 7 +DEAL::Copy particle location: 0.300000 0.500000 +DEAL::Copy particle reference location: 0.200000 0.400000 +DEAL::Copy particle index: 7 +DEAL::Moved particle location: 0.300000 0.500000 +DEAL::Moved particle reference location: 0.200000 0.400000 +DEAL::Moved particle index: 7 +DEAL::Original particle location: 0.300000 0.500000 +DEAL::Original particle reference location: 0.200000 0.400000 +DEAL::Original particle index: 7 +DEAL::OK diff --git a/tests/particles/particle_03.cc b/tests/particles/particle_03.cc new file mode 100644 index 0000000000..06d2685f85 --- /dev/null +++ b/tests/particles/particle_03.cc @@ -0,0 +1,75 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// Like particle_02, but with particle properties. + +#include "../tests.h" +#include +#include +#include + + +template +void test () +{ + { + const unsigned int n_properties_per_particle = 3; + Particles::PropertyPool pool(n_properties_per_particle); + + Point<2> position; + position(0) = 0.3; + position(1) = 0.5; + + Point<2> reference_position; + reference_position(0) = 0.2; + reference_position(1) = 0.4; + + const Particles::types::particle_index index(7); + + std::vector properties = {0.15,0.45,0.75}; + + Particles::Particle<2> particle(position,reference_position,index); + particle.set_property_pool(pool); + particle.set_properties(properties); + + deallog << "Particle properties: " + << std::vector(particle.get_properties().begin(),particle.get_properties().end()) + << std::endl; + + const Particles::Particle<2> copy(particle); + + deallog << "Copy particle properties: " + << std::vector(copy.get_properties().begin(),copy.get_properties().end()) + << std::endl; + + const Particles::Particle<2> moved_particle(std::move(particle)); + + deallog << "Moved particle properties: " + << std::vector(moved_particle.get_properties().begin(),moved_particle.get_properties().end()) + << std::endl; + } + + deallog << "OK" << std::endl; +} + + + +int main () +{ + initlog(); + test<2>(); +} diff --git a/tests/particles/particle_03.output b/tests/particles/particle_03.output new file mode 100644 index 0000000000..558d9d51af --- /dev/null +++ b/tests/particles/particle_03.output @@ -0,0 +1,5 @@ + +DEAL::Particle properties: 0.150000 0.450000 0.750000 +DEAL::Copy particle properties: 0.150000 0.450000 0.750000 +DEAL::Moved particle properties: 0.150000 0.450000 0.750000 +DEAL::OK