Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
particle_coupling.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2024 by Luca Heltai
4//
5// This file is part of the reduced_lagrange_multipliers application, based on
6// the deal.II library.
7//
8// The reduced_lagrange_multipliers application is free software; you can use
9// it, redistribute it, and/or modify it under the terms of the Apache-2.0
10// License WITH LLVM-exception as published by the Free Software Foundation;
11// either version 3.0 of the License, or (at your option) any later version. The
12// full text of the license can be found in the file LICENSE.md at the top level
13// of the reduced_lagrange_multipliers distribution.
14//
15// ---------------------------------------------------------------------
16
17#include "particle_coupling.h"
18template <int dim>
24
25template <int dim>
31
32
33
34template <int dim>
35void
36ParticleCoupling<dim>::output_particles(const std::string &output_name) const
37{
38 Particles::DataOut<dim> particles_out;
39 particles_out.build_patches(particles);
40 particles_out.write_vtu_in_parallel(output_name, mpi_communicator);
41}
42
43
44
45template <int dim>
46void
49 const Mapping<dim> &mapp)
50{
51 tria_background = &tria;
52 mapping = &mapp;
53 particles.initialize(*tria_background, *mapping, 1);
54 mpi_communicator = tria_background->get_communicator();
56 {
57 std::vector<BoundingBox<dim>> all_boxes;
58 all_boxes.reserve(tria.n_locally_owned_active_cells());
59 for (const auto &cell : tria.active_cell_iterators())
60 if (cell->is_locally_owned())
61 all_boxes.emplace_back(cell->bounding_box());
62 const auto tree = pack_rtree(all_boxes);
63 const auto local_boxes =
64 extract_rtree_level(tree, par.rtree_extraction_level);
65
66 global_bounding_boxes =
67 Utilities::MPI::all_gather(mpi_communicator, local_boxes);
68 }
69
70 Assert(!global_bounding_boxes.empty(),
71 ExcInternalError(
72 "I was expecting the "
73 "global_bounding_boxes to be filled at this stage. "
74 "Make sure you fill this vector before trying to use it "
75 "here. Bailing out."));
76}
77
78
79
80template <int dim>
81std::vector<std::vector<BoundingBox<dim>>>
87
88
89template <int dim>
92{
93 return particles;
94}
95
96
97
98template <int dim>
99std::map<unsigned int, IndexSet>
101 const std::vector<Point<dim>> &points,
102 const std::vector<std::vector<double>> &properties)
103{
104 AssertThrow(tria_background, ExcNotInitialized());
105 auto local_indices_map =
106 particles.insert_global_particles(points,
108 properties);
109 return local_indices_map;
111
112// Explicit instantiations
113
114template class ParticleCouplingParameters<2>;
115template class ParticleCouplingParameters<3>;
117template class ParticleCoupling<2>;
118template class ParticleCoupling<3>;
void write_vtu_in_parallel(const std::string &filename, const MPI_Comm comm) const
ParameterAcceptor(const std::string &section_name="")
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", ParameterHandler &prm_=prm, const Patterns::PatternBase &pattern=*Patterns::Tools::Convert< ParameterType >::to_pattern())
void build_patches(const Particles::ParticleHandler< dim, spacedim > &particles, const std::vector< std::string > &data_component_names={}, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretations={})
SmartPointer< const Mapping< dim > > mapping
Smart pointer to the mapping associated with the triangulation.
std::vector< std::vector< BoundingBox< dim > > > global_bounding_boxes
void output_particles(const std::string &output_name) const
Outputs the current state of the particles to a file.
std::map< unsigned int, IndexSet > insert_points(const std::vector< Point< dim > > &points, const std::vector< std::vector< double > > &properties={})
std::vector< std::vector< BoundingBox< dim > > > get_global_bounding_boxes() const
const Particles::ParticleHandler< dim > & get_particles() const
MPI_Comm mpi_communicator
Get the MPI communicator associated with the triangulation.
SmartPointer< const parallel::TriangulationBase< dim > > tria_background
Smart pointer to the background triangulation.
Particles::ParticleHandler< dim > particles
Handler for managing particles in the simulation.
const ParticleCouplingParameters< dim > & par
Parameters for particle coupling.
void initialize_particle_handler(const parallel::TriangulationBase< dim > &tria_background, const Mapping< dim > &mapping=StaticMappingQ1< dim >::mapping)
ParticleCoupling(const ParticleCouplingParameters< dim > &par)
Constructor.
Stores parameters related to particle coupling in a simulation.
ParticleCouplingParameters()
Constructor that initializes the parameters.
unsigned int n_locally_owned_active_cells() const
IteratorRange< active_cell_iterator > active_cell_iterators() const
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
std::vector< T > all_gather(const MPI_Comm comm, const T &object_to_send)
std::vector< BoundingBox< boost::geometry::dimension< typename Rtree::indexable_type >::value > > extract_rtree_level(const Rtree &tree, const unsigned int level)
RTree< typename LeafTypeIterator::value_type, IndexType, IndexableGetter > pack_rtree(const LeafTypeIterator &begin, const LeafTypeIterator &end)