Reduced Lagrange Multipliers
 
Loading...
Searching...
No Matches
immersed_repartitioner.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
18
20
21template <int dim, int spacedim>
26
27template <int dim, int spacedim>
30 const Triangulation<dim, spacedim> &tria_immersed) const
31{
32 // 1) collect centers of immeresed mesh
33 std::vector<Point<spacedim>> points;
34
35 for (const auto &cell : tria_immersed.active_cell_iterators())
36 if (cell->is_locally_owned())
37 points.push_back(cell->center());
38
39 // 2) determine owner on background mesh
41 Vector<double> ranks(tria_background.n_active_cells());
42 ranks = Utilities::MPI::this_mpi_process(tria_background.get_communicator());
43
44 const auto point_ranks =
47 ranks,
48 points,
49 rpe,
51 0);
52
53 const auto tria =
55 &tria_immersed);
56
57 Assert(tria, ExcNotImplemented());
58
59 // 3) set partitioning
61 tria->global_active_cell_index_partitioner().lock());
62
63 unsigned int counter = 0;
64 for (const auto &cell : tria_immersed.active_cell_iterators())
65 if (cell->is_locally_owned())
66 partition[cell->global_active_cell_index()] = point_ranks[counter++];
67
68 partition.update_ghost_values();
70 return partition;
71}
72
73
74template class ImmersedRepartitioner<1, 1>;
75template class ImmersedRepartitioner<1, 2>;
76template class ImmersedRepartitioner<1, 3>;
77template class ImmersedRepartitioner<2, 2>;
78template class ImmersedRepartitioner<2, 3>;
79template class ImmersedRepartitioner<3, 3>;
A class implementing a repartitioning policy for immersed triangulations.
ImmersedRepartitioner(const Triangulation< spacedim > &tria_background)
Constructor for the ImmersedRepartitioner class.
const MappingQ1< spacedim > mapping
A mapping object for the background triangulation.
virtual LinearAlgebra::distributed::Vector< double > partition(const Triangulation< dim, spacedim > &tria_immersed) const override
Repartition the given immersed triangulation.
const Triangulation< spacedim > & tria_background
A reference to the background triangulation used for repartitioning.
IteratorRange< active_cell_iterator > active_cell_iterators() const
#define Assert(cond, exc)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
std::vector< typename FEPointEvaluation< n_components, dim, spacedim, typename VectorType::value_type >::value_type > point_values(const Mapping< dim > &mapping, const MeshType< dim, spacedim > &mesh, const VectorType &vector, const std::vector< Point< spacedim > > &evaluation_points, Utilities::MPI::RemotePointEvaluation< dim, spacedim > &cache, const EvaluationFlags::EvaluationFlags flags=EvaluationFlags::avg, const unsigned int first_selected_component=0)