File: assemble_vector.h

package info (click to toggle)
dolfinx-mpc 0.9.3-1
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 1,188 kB
  • sloc: python: 7,263; cpp: 5,462; makefile: 69; sh: 4
file content (115 lines) | stat: -rw-r--r-- 4,414 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
// Copyright (C) 2021 Jorgen S. Dokken, Nathan Sime, and Connor D. Pierce
//
// This file is part of DOLFINX_MPC
//
// SPDX-License-Identifier:    MIT

#pragma once

#include "MultiPointConstraint.h"
#include "assemble_utils.h"
#include <dolfinx/fem/DirichletBC.h>
#include <dolfinx/fem/Form.h>
#include <functional>
using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
    const std::int32_t,
    MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;

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

/// Given a local element vector, move all slave contributions to the global
/// (local to process) vector.
/// @param [in, out] b The global (local to process) vector
/// @param [in, out] b_local The local element vector
/// @param [in] b_local_copy Copy of the local element vector
/// @param [in] dofs The cell dofs (blocked)
/// @param [in] num_dofs The number of degrees of freedom in the local vector
/// @param [in] bs The element block size
/// @param [in] is_slave Vector indicating if a dof is a slave
/// @param [in] slaves The slave dofs (local to process)
/// @param [in] masters Adjacency list with master dofs
/// @param [in] coeffs Adjacency list with the master coefficients
template <typename T>
void modify_mpc_vec(
    const std::span<T>& b, const std::span<T>& b_local,
    const std::span<T>& b_local_copy, std::span<const std::int32_t> dofs,
    const int num_dofs, const int bs, std::span<const std::int8_t> is_slave,
    std::span<const std::int32_t> slaves,
    const std::shared_ptr<const dolfinx::graph::AdjacencyList<std::int32_t>>&
        masters,
    const std::shared_ptr<const dolfinx::graph::AdjacencyList<T>>& coeffs)
{

  // NOTE: Should this be moved into the MPC constructor?
  // Get local index of slaves in cell
  std::vector<std::int32_t> local_index
      = compute_local_slave_index(slaves, num_dofs, bs, dofs, is_slave);

  // Move contribution from each slave to corresponding master dof
  for (std::size_t i = 0; i < local_index.size(); i++)
  {
    auto masters_i = masters->links(slaves[i]);
    auto coeffs_i = coeffs->links(slaves[i]);
    assert(masters_i.size() == coeffs_i.size());
    for (std::size_t j = 0; j < masters_i.size(); j++)
    {
      if constexpr (std::is_scalar_v<T>)
        // Use standard transpose for type double
        b[masters_i[j]] += coeffs_i[j] * b_local_copy[local_index[i]];
      else
        // Use Hermitian transpose for type std::complex<double>
        b[masters_i[j]]
            += std::conj(coeffs_i[j]) * b_local_copy[local_index[i]];
      b_local[local_index[i]] = 0;
    }
  }
}

/// Assemble a linear form into a vector
/// @param[in] b The vector to be assembled. It will not be zeroed before
/// assembly.
/// @param[in] L The linear forms to assemble into b
/// @param[in] mpc The multi-point constraint
void assemble_vector(
    std::span<double> b, const dolfinx::fem::Form<double>& L,
    const std::shared_ptr<
        const dolfinx_mpc::MultiPointConstraint<double, double>>& mpc);

/// Assemble a linear form into a vector
/// @param[in] b The vector to be assembled. It will not be zeroed before
/// assembly.
/// @param[in] L The linear forms to assemble into b
/// @param[in] mpc The multi-point constraint
void assemble_vector(
    std::span<std::complex<double>> b,
    const dolfinx::fem::Form<std::complex<double>>& L,
    const std::shared_ptr<
        const dolfinx_mpc::MultiPointConstraint<std::complex<double>, double>>&
        mpc);

/// Assemble a linear form into a vector
/// @param[in] b The vector to be assembled. It will not be zeroed before
/// assembly.
/// @param[in] L The linear forms to assemble into b
/// @param[in] mpc The multi-point constraint
void assemble_vector(
    std::span<float> b, const dolfinx::fem::Form<float>& L,
    const std::shared_ptr<
        const dolfinx_mpc::MultiPointConstraint<float, float>>& mpc);

/// Assemble a linear form into a vector
/// @param[in] b The vector to be assembled. It will not be zeroed before
/// assembly.
/// @param[in] L The linear forms to assemble into b
/// @param[in] mpc The multi-point constraint
void assemble_vector(
    std::span<std::complex<float>> b,
    const dolfinx::fem::Form<std::complex<float>>& L,
    const std::shared_ptr<
        const dolfinx_mpc::MultiPointConstraint<std::complex<float>, float>>&
        mpc);

} // namespace dolfinx_mpc