File: MultiPointConstraint.h

package info (click to toggle)
dolfinx-mpc 0.9.3-1
  • links: PTS, VCS
  • area: main
  • in suites:
  • size: 1,188 kB
  • sloc: python: 7,263; cpp: 5,462; makefile: 69; sh: 4
file content (223 lines) | stat: -rw-r--r-- 8,202 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
// Copyright (C) 2019-2021 Jorgen S. Dokken
//
// This file is part of DOLFINX_MPC
//
// SPDX-License-Identifier:    MIT

#pragma once

#include "mpc_helpers.h"
#include <dolfinx/common/IndexMap.h>
#include <dolfinx/common/Timer.h>
#include <dolfinx/common/log.h>
#include <dolfinx/fem/DofMap.h>
#include <dolfinx/fem/FunctionSpace.h>
#include <dolfinx/graph/AdjacencyList.h>
#include <dolfinx/mesh/Mesh.h>
#include <iostream>

namespace dolfinx_mpc
{
template <typename T, std::floating_point U>
class MultiPointConstraint

{

public:
  /// Create contact constraint
  ///
  /// @param[in] V The function space
  /// @param[in] slaves List of local slave dofs
  /// @param[in] masters Array of all masters
  /// @param[in] coeffs Coefficients corresponding to each master
  /// @param[in] owners Owners for each master
  /// @param[in] offsets Offsets for masters
  /// @tparam The floating type of the mesh
  MultiPointConstraint(std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> V,
                       std::span<const std::int32_t> slaves,
                       std::span<const std::int64_t> masters,
                       std::span<const T> coeffs,
                       std::span<const std::int32_t> owners,
                       std::span<const std::int32_t> offsets)
      : _slaves(), _is_slave(), _cell_to_slaves_map(), _num_local_slaves(),
        _master_map(), _coeff_map(), _owner_map(), _mpc_constants(), _V()
  {
    assert(slaves.size() == offsets.size() - 1);
    assert(masters.size() == coeffs.size());
    assert(coeffs.size() == owners.size());
    assert(offsets.back() == owners.size());

    // Create list indicating which dofs on the process are slaves
    const dolfinx::fem::DofMap& dofmap = *(V->dofmap());
    const std::int32_t num_dofs_local
        = dofmap.index_map_bs()
          * (dofmap.index_map->size_local() + dofmap.index_map->num_ghosts());
    std::vector<std::int8_t> _slave_data(num_dofs_local, 0);
    _mpc_constants = std::vector<T>(num_dofs_local, 0);
    for (auto dof : slaves)
    {
      _slave_data[dof] = 1;
      // FIXME: Add input vector for this data
      _mpc_constants[dof] = 1;
    }
    _is_slave = std::move(_slave_data);

    // Create a map for cells owned by the process to the slaves
    _cell_to_slaves_map = create_cell_to_dofs_map(*V, slaves);

    // Create adjacency list with all local dofs, where the slave dofs maps to
    // its masters
    std::vector<std::int32_t> _num_masters(num_dofs_local);
    std::ranges::fill(_num_masters, 0);
    for (std::int32_t i = 0; i < slaves.size(); i++)
      _num_masters[slaves[i]] = offsets[i + 1] - offsets[i];
    std::vector<std::int32_t> masters_offsets(num_dofs_local + 1);
    masters_offsets[0] = 0;
    std::inclusive_scan(_num_masters.begin(), _num_masters.end(),
                        masters_offsets.begin() + 1);

    // Reuse num masters as fill position array
    std::ranges::fill(_num_masters, 0);
    std::vector<std::int64_t> _master_data(masters.size());
    std::vector<T> _coeff_data(masters.size());
    std::vector<std::int32_t> _owner_data(masters.size());
    /// Create adjacency lists spanning all local dofs mapping to master dofs,
    /// its owner and the corresponding coefficient
    for (std::size_t i = 0; i < slaves.size(); i++)
    {
      for (std::int32_t j = 0; j < offsets[i + 1] - offsets[i]; j++)
      {
        _master_data[masters_offsets[slaves[i]] + _num_masters[slaves[i]]]
            = masters[offsets[i] + j];
        _coeff_data[masters_offsets[slaves[i]] + _num_masters[slaves[i]]]
            = coeffs[offsets[i] + j];
        _owner_data[masters_offsets[slaves[i]] + _num_masters[slaves[i]]]
            = owners[offsets[i] + j];
        _num_masters[slaves[i]]++;
      }
    }
    _coeff_map = std::make_shared<dolfinx::graph::AdjacencyList<T>>(
        _coeff_data, masters_offsets);
    _owner_map = std::make_shared<dolfinx::graph::AdjacencyList<std::int32_t>>(
        _owner_data, masters_offsets);

    // Create a vector containing all the slave dofs (sorted)
    std::vector<std::int32_t> sorted_slaves(slaves.size());
    std::int32_t c = 0;
    for (std::size_t i = 0; i < _is_slave.size(); i++)
      if (_is_slave[i])
        sorted_slaves[c++] = i;
    _slaves = std::move(sorted_slaves);

    const std::int32_t num_local
        = dofmap.index_map_bs() * dofmap.index_map->size_local();
    auto it = std::ranges::lower_bound(_slaves, num_local);
    _num_local_slaves = std::distance(_slaves.begin(), it);

    // Create new function space with extended index map
    _V = std::make_shared<const dolfinx::fem::FunctionSpace<U>>(
        create_extended_functionspace(*V, _master_data, _owner_data));

    // Map global masters to local index in extended function space
    std::vector<std::int32_t> masters_local
        = map_dofs_global_to_local<U>(*_V, _master_data);
    _master_map = std::make_shared<dolfinx::graph::AdjacencyList<std::int32_t>>(
        masters_local, masters_offsets);
  }
  //-----------------------------------------------------------------------------
  /// Backsubstitute slave/master constraint for a given function
  void backsubstitution(std::span<T> vector)
  {
    for (auto slave : _slaves)
    {
      // Zero out intial data in slave dof
      vector[slave] = 0.0;
      // Accumulate master contributions
      auto masters = _master_map->links(slave);
      auto coeffs = _coeff_map->links(slave);
      assert(masters.size() == coeffs.size());
      for (std::size_t k = 0; k < masters.size(); ++k)
        vector[slave]
            += coeffs[k] * vector[masters[k]]; //+ _mpc_constants[slave];
    }
  };

  /// Homogenize slave DoFs (particularly useful for nonlinear problems)
  void homogenize(std::span<T> vector) const
  {
    for (auto slave : _slaves)
      vector[slave] = 0.0;
  };

  /// Return map from cell to slaves contained in that cell
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>
  cell_to_slaves() const
  {
    return _cell_to_slaves_map;
  }
  /// Return map from slave to masters (local_index)
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>
  masters() const
  {
    return _master_map;
  }

  /// Return map from slave to coefficients
  std::shared_ptr<const dolfinx::graph::AdjacencyList<T>> coefficients() const
  {
    return _coeff_map;
  }

  /// Return map from slave to masters (global index)
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>
  owners() const
  {
    return _owner_map;
  }

  /// Return of local dofs + num ghosts indicating if a dof is a slave
  std::span<const std::int8_t> is_slave() const
  {
    return std::span<const std::int8_t>(_is_slave);
  }

  /// Return the constant values for the constraint
  const std::vector<T>& constant_values() const { return _mpc_constants; }

  /// Return an array of all slave indices (sorted and local to process)
  const std::vector<std::int32_t>& slaves() const { return _slaves; }

  /// Return number of slaves owned by process
  const std::int32_t num_local_slaves() const { return _num_local_slaves; }

  /// Return the MPC FunctionSpace
  std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> function_space() const
  {
    return _V;
  }

private:
  // MPC function space
  std::shared_ptr<const dolfinx::fem::FunctionSpace<U>> _V;

  // Array including all slaves (local + ghosts)
  std::vector<std::int32_t> _slaves;
  std::vector<std::int8_t> _is_slave;

  std::vector<T> _mpc_constants;

  // Map from slave cell to index in _slaves for a given slave cell
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>
      _cell_to_slaves_map;

  // Number of slaves owned by the process
  std::int32_t _num_local_slaves;
  // Map from slave (local to process) to masters (local to process)
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>
      _master_map;
  // Map from slave (local to process)to coefficients
  std::shared_ptr<const dolfinx::graph::AdjacencyList<T>> _coeff_map;
  // Map from slave( local to process) to rank of process owning master
  std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>> _owner_map;
};
} // namespace dolfinx_mpc