GeneralizedTAOSolver.h 12.3 KB
Newer Older
webmanue's avatar
webmanue committed
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
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
// © 2020, 2022 ETH Zurich, Mechanics and Materials Lab
// © 2020 California Institute of Technology
//
// Licensed under the Apache License, Version 2.0 (the "License");
// you may not use this file except in compliance with the License.
// You may obtain a copy of the License at
//
// http://www.apache.org/licenses/LICENSE-2.0
//
// Unless required by applicable law or agreed to in writing, software
// distributed under the License is distributed on an "AS IS" BASIS,
// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
// See the License for the specific language governing permissions and
// limitations under the License.

#pragma once

#ifndef AE108_PETSC_COMPLEX

#include "ae108/assembly/AssemblerTypeTraits.h"
#include "ae108/cpppetsc/MeshBoundaryCondition.h"
#include "ae108/cpppetsc/TAOSolver.h"
#include "ae108/cpppetsc/TaggedVector.h"
#include <cassert>
#include <functional>
#include <mpi.h>
#include <petscsys.h>
#include <vector>

namespace ae108 {
namespace solve {

template <class Assembler> class GeneralizedTAOSolver {
public:
  using policy_type = typename assembly::PolicyTypeTrait<Assembler>::type;
  using mesh_type = typename assembly::MeshTypeTrait<Assembler>::type;
  using size_type = typename mesh_type::size_type;
  using value_type = typename mesh_type::value_type;
  using real_type = typename mesh_type::real_type;
  using vector_type = typename mesh_type::vector_type;
  using matrix_type = typename mesh_type::matrix_type;

  /**
   * @param mesh A valid nonzero pointer to a mesh_type instance.
   */
  explicit GeneralizedTAOSolver(const mesh_type *mesh);

  struct BoundaryConditionContainer {
    /**
     * @brief The dimension of the output of `residual`.
     */
    size_type size;

    /**
     * @brief The solver attempts to reach a residual of zero.
     */
    typename cpppetsc::TAOSolver<policy_type>::GradientFunctor residual;

    /**
     * @brief The jacobian of the residual.
     */
    typename cpppetsc::TAOSolver<policy_type>::HessianFunctor jacobian;
  };

  /**
   * @brief Calls computeSolution with the default assembler calls (ie. calls
   * assembleEnergy, assembleForceVector or assembleStiffnessMatrix). Passes on
   * the rest of the arguments.
   *
   * @param boundaryConditions The nonlinear boundary conditions to apply.
   * @param initialGuess The global vector to start iterating at.
   * @param time This will be used to configure the assembler.
   * @param assembler Valid nonzero pointer.
   */
  cpppetsc::distributed<vector_type>
  computeSolution(const BoundaryConditionContainer &boundaryConditions,
                  cpppetsc::distributed<vector_type> initialGuess,
                  const double time, const Assembler *const assembler) const;

  using LocalEnergyAssembler = std::function<void(
      const cpppetsc::local<vector_type> &, double, real_type *)>;

  using LocalForceVectorAssembler =
      std::function<void(const cpppetsc::local<vector_type> &, double,
                         cpppetsc::local<vector_type> *)>;

  using LocalStiffnessMatrixAssembler = std::function<void(
      const cpppetsc::local<vector_type> &, double, matrix_type *)>;

  /**
   * @brief Calls computeSolution wrapping up the local assembler calls into
   * distributed assembler calls. Passes on the rest of the arguments.
   *
   * @param boundaryConditions The nonlinear boundary conditions to apply.
   * Note that empty functions (residual, jacobian) are interpreted as functions
   * that do not change the output.
   * @param initialGuess The global vector to start iterating at.
   * @param time This will be used to configure the assembler.
   * @param assembleEnergy A valid callable. It will be called to assemble
   * the local energy.
   * @param assembleForceVector A valid callable. It will be called to assemble
   * the local force vector.
   * @param assembleStiffnessMatrix A valid callable. It will be called to
   * assemble the stiffness matrix.
   */
  cpppetsc::distributed<vector_type>
  computeSolution(const BoundaryConditionContainer &boundaryConditions,
                  cpppetsc::distributed<vector_type> initialGuess,
                  const double time, LocalEnergyAssembler assembleEnergy,
                  LocalForceVectorAssembler assembleForceVector,
                  LocalStiffnessMatrixAssembler assemblerStiffnessMatrix) const;

  using DistributedEnergyAssembler = std::function<void(
      const cpppetsc::distributed<vector_type> &, double, real_type *)>;

  using DistributedForceVectorAssembler =
      std::function<void(const cpppetsc::distributed<vector_type> &, double,
                         cpppetsc::distributed<vector_type> *)>;

  using DistributedStiffnessMatrixAssembler = std::function<void(
      const cpppetsc::distributed<vector_type> &, double, matrix_type *)>;

  /**
   * @brief The goal of this function is to solve E(x) = min under the condition
   * that the residual is equal to zero.
   *
   * @remark Does not update internal variables.
   *
   * @param boundaryConditions The nonlinear boundary conditions to apply.
   * @param initialGuess The global vector to start iterating at.
   * @param time This will be used to configure the assembler.
   * @param assembleEnergy A valid callable. It will be called to assemble
   * the distributed energy.
   * @param assembleForceVector A valid callable. It will be called to assemble
   * the distributed force vector.
   * @param assembleStiffnessMatrix A valid callable. It will be called to
   * assemble the stiffness matrix.
   */
  cpppetsc::distributed<vector_type> computeSolution(
      const BoundaryConditionContainer &boundaryConditions,
      cpppetsc::distributed<vector_type> initialGuess, const double time,
      DistributedEnergyAssembler assembleEnergy,
      DistributedForceVectorAssembler assembleForceVector,
      DistributedStiffnessMatrixAssembler assemblerStiffnessMatrix) const;

private:
  const mesh_type *_mesh;
}; // namespace solve
} // namespace solve
} // namespace ae108

/********************************************************************
 *  implementations
 *******************************************************************/

namespace ae108 {
namespace solve {
template <class Assembler>
GeneralizedTAOSolver<Assembler>::GeneralizedTAOSolver(const mesh_type *mesh)
    : _mesh(mesh) {}

template <class Assembler>
cpppetsc::distributed<typename GeneralizedTAOSolver<Assembler>::vector_type>
GeneralizedTAOSolver<Assembler>::computeSolution(
    const BoundaryConditionContainer &boundaryConditions,
    cpppetsc::distributed<vector_type> initialGuess, const double time,
    const Assembler *const assembler) const {
  assert(assembler);
  return computeSolution(
      boundaryConditions, std::move(initialGuess), time,
      LocalEnergyAssembler(
          [assembler](const cpppetsc::local<vector_type> &localDisplacements,
                      const double time, real_type *const localEnergy) {
            assembler->assembleEnergy(localDisplacements, time, localEnergy);
          }),
      LocalForceVectorAssembler(
          [assembler](const cpppetsc::local<vector_type> &localDisplacements,
                      const double time,
                      cpppetsc::local<vector_type> *const localForces) {
            assembler->assembleForceVector(localDisplacements, time,
                                           localForces);
          }),
      LocalStiffnessMatrixAssembler(
          [assembler](const cpppetsc::local<vector_type> &localDisplacements,
                      const double time, matrix_type *const output) {
            assembler->assembleStiffnessMatrix(localDisplacements, time,
                                               output);
          }));
}

template <class Assembler>
cpppetsc::distributed<typename GeneralizedTAOSolver<Assembler>::vector_type>
GeneralizedTAOSolver<Assembler>::computeSolution(
    const BoundaryConditionContainer &boundaryConditions,
    cpppetsc::distributed<vector_type> initialGuess, const double time,
    LocalEnergyAssembler assembleEnergy,
    LocalForceVectorAssembler assembleForceVector,
    LocalStiffnessMatrixAssembler assembleStiffnessMatrix) const {
  assert(assembleEnergy);
  assert(assembleForceVector);
  assert(assembleStiffnessMatrix);

  const auto &mesh = *_mesh;

  auto localDisplacements = vector_type::fromLocalMesh(mesh);
  auto localForces = vector_type::fromLocalMesh(mesh);

  const auto distributedEnergyAssembler =
      [&assembleEnergy, &mesh,
       &localDisplacements](const cpppetsc::distributed<vector_type> &input,
                            const double time, real_type *const output) {
        mesh.copyToLocalVector(input, &localDisplacements);

        *output = real_type{};
        assembleEnergy(localDisplacements, time, output);

        policy_type::handleError(MPI_Allreduce(MPI_IN_PLACE, output, 1,
                                               MPIU_REAL, MPIU_SUM,
                                               policy_type::communicator()));
      };

  const auto distributedForceVectorAssembler =
      [&assembleForceVector, &mesh, &localDisplacements, &localForces](
          const cpppetsc::distributed<vector_type> &input, const double time,
          cpppetsc::distributed<vector_type> *const output) {
        mesh.copyToLocalVector(input, &localDisplacements);

        localForces.unwrap().setZero();
        assembleForceVector(localDisplacements, time, &localForces);

        mesh.addToGlobalVector(localForces, output);
      };

  const auto distributedStiffnessMatrixAssembler =
      [&assembleStiffnessMatrix, &mesh,
       &localDisplacements](const cpppetsc::distributed<vector_type> &input,
                            const double time, matrix_type *const output) {
        mesh.copyToLocalVector(input, &localDisplacements);

        assembleStiffnessMatrix(localDisplacements, time, output);
        output->finalize();
      };

  return computeSolution(boundaryConditions, std::move(initialGuess), time,
                         distributedEnergyAssembler,
                         distributedForceVectorAssembler,
                         distributedStiffnessMatrixAssembler);
}

template <class Assembler>
cpppetsc::distributed<typename GeneralizedTAOSolver<Assembler>::vector_type>
GeneralizedTAOSolver<Assembler>::computeSolution(
    const BoundaryConditionContainer &boundaryConditions,
    cpppetsc::distributed<vector_type> initialGuess, const double time,
    DistributedEnergyAssembler assembleEnergy,
    DistributedForceVectorAssembler assembleForceVector,
    DistributedStiffnessMatrixAssembler assembleStiffnessMatrix) const {
  assert(assembleEnergy);
  assert(assembleForceVector);
  assert(assembleStiffnessMatrix);

  using Solver = cpppetsc::TAOSolver<policy_type>;
  auto matrix = matrix_type::fromMesh(*_mesh);

  const auto global = matrix.size();
  const auto local = matrix.localSize();

webmanue's avatar
webmanue committed
268
  Solver solver{std::move(matrix), Solver::Type::almm};
webmanue's avatar
webmanue committed
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311

  solver.setConstraints(typename Solver::EqualityConstraints{
      boundaryConditions.residual
          ? boundaryConditions.residual
          : [](const cpppetsc::distributed<vector_type> &,
               cpppetsc::distributed<vector_type> *) {},
      boundaryConditions.jacobian
          ? boundaryConditions.jacobian
          : [](const cpppetsc::distributed<vector_type> &, matrix_type *) {},
      {
          cpppetsc::distributed<vector_type>(boundaryConditions.size),
          matrix_type(typename matrix_type::LocalRows(PETSC_DECIDE),
                      typename matrix_type::LocalCols(local.second),
                      typename matrix_type::GlobalRows(boundaryConditions.size),
                      typename matrix_type::GlobalCols(global.second)),
      },
  });

  return solver.solve(
      [time, &assembleEnergy](const cpppetsc::distributed<vector_type> &input,
                              double *const output) {
        *output = real_type{};
        assembleEnergy(input, time, output);
      },
      [time,
       &assembleForceVector](const cpppetsc::distributed<vector_type> &input,
                             cpppetsc::distributed<vector_type> *const output) {
        output->unwrap().setZero();
        assembleForceVector(input, time, output);
      },
      [time, &assembleStiffnessMatrix](
          const cpppetsc::distributed<vector_type> &input,
          matrix_type *const output) {
        output->setZero();
        assembleStiffnessMatrix(input, time, output);
      },
      std::move(initialGuess));
}

} // namespace solve
} // namespace ae108

#endif