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
268
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
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457 | // *****************************************************************************
/*!
\file src/Inciter/ChoCG.hpp
\copyright 2012-2015 J. Bakosi,
2016-2018 Los Alamos National Security, LLC.,
2019-2021 Triad National Security, LLC.,
2022-2025 J. Bakosi
All rights reserved. See the LICENSE file for details.
\brief ChoCG: Projection-based solver for incompressible flow
*/
// *****************************************************************************
#pragma once
#include <vector><--- Include file: not found. Please note: Cppcheck does not need standard library headers to get proper results.
#include <map><--- Include file:
#include "Types.hpp"
#include "Fields.hpp"
#include "Table.hpp"
#include "DerivedData.hpp"<--- Include file: "DerivedData.hpp" not found.
#include "NodeDiagnostics.hpp"
#include "PUPUtil.hpp"
#include "ConjugateGradients.hpp"<--- Include file: "ConjugateGradients.hpp" not found.
#include "NoWarning/chocg.decl.h"<--- Include file: "NoWarning/chocg.decl.h" not found.
namespace inciter {
//! ChoCG Charm++ chare array used to advance PDEs in time with ChoCG
class ChoCG : public CBase_ChoCG {
public:
#if defined(__clang__)
#pragma clang diagnostic push
#pragma clang diagnostic ignored "-Wunused-parameter"
#pragma clang diagnostic ignored "-Wdeprecated-declarations"
#elif defined(STRICT_GNUC)
#pragma GCC diagnostic push
#pragma GCC diagnostic ignored "-Wunused-parameter"
#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
#elif defined(__INTEL_COMPILER)
#pragma warning( push )
#pragma warning( disable: 1478 )
#endif
// Include Charm++ SDAG code. See http://charm.cs.illinois.edu/manuals/html/
// charm++/manual.html, Sec. "Structured Control Flow: Structured Dagger".
ChoCG_SDAG_CODE
#if defined(__clang__)
#pragma clang diagnostic pop
#elif defined(STRICT_GNUC)
#pragma GCC diagnostic pop
#elif defined(__INTEL_COMPILER)
#pragma warning( pop )
#endif
//! Constructor
explicit ChoCG( const CProxy_Discretization& disc,
const tk::CProxy_ConjugateGradients& cgpre,
const tk::CProxy_ConjugateGradients& cgmom,
const std::map< int, std::vector< std::size_t > >& bface,
const std::map< int, std::vector< std::size_t > >& bnode,
const std::vector< std::size_t >& triinpoel );
#if defined(__clang__)
#pragma clang diagnostic push
#pragma clang diagnostic ignored "-Wundefined-func-template"
#endif
//! Migrate constructor
// cppcheck-suppress uninitMemberVar
explicit ChoCG( CkMigrateMessage* m ) : CBase_ChoCG( m ) {}
#if defined(__clang__)
#pragma clang diagnostic pop
#endif
//! Configure Charm++ custom reduction types initiated from this chare array
static void registerReducers();
//! Return from migration
void ResumeFromSync() override;
//! Start setup for solution
void setup( tk::real v );
//! Initialize Poisson solve
void pinit();
//! Solve Poisson equation
void psolve();
//! Continue after Poisson solve
void psolved();
//! Solve momentum/transport equations
void msolve();
//! Continue after momentum/transport ssolve
void msolved();
// Start time stepping
void start();
//! Advance equations to next time step
void advance( tk::real newdt );
//! Compute righ-hand side vector of transport equations
void rhs();
//! Evaluate diagnostics
void diag();
//! Start (re-)computing domain and boundary integrals
void feop();
//! Receive contributions to boundary point normals on chare-boundaries
void comnorm( const std::unordered_map< int,
std::unordered_map< std::size_t, std::array<tk::real,4> > >& inbnd );
//! Receive contributions to velocity gradients
void comvgrad( const std::unordered_map< std::size_t,
std::vector< tk::real > >& ingrad );
//! Receive contributions to momentum flux on chare-boundaries
void comflux( const std::unordered_map< std::size_t,
std::vector< tk::real > >& influx );
//! Receive contributions to conjugate gradients solution gradient
void comsgrad( const std::unordered_map< std::size_t,
std::vector< tk::real > >& ingrad );
//! Receive contributions to pressure gradient
void compgrad( const std::unordered_map< std::size_t,
std::vector< tk::real > >& ingrad );
//! Receive contributions to right-hand side vector on chare-boundaries
void comrhs( const std::unordered_map< std::size_t,
std::vector< tk::real > >& inrhs );
//! Receive contributions to velocity divergence on chare-boundaries
void comdiv( const std::unordered_map< std::size_t, tk::real >& indiv );
//! Evaluate residuals
void evalres( const std::vector< tk::real >& l2res );
//! Receive new mesh from Refiner
void resizePostAMR(
const std::vector< std::size_t >& ginpoel,
const tk::UnsMesh::Chunk& chunk,
const tk::UnsMesh::Coords& coord,
const std::unordered_map< std::size_t, tk::UnsMesh::Edge >& addedNodes,
const std::unordered_map< std::size_t, std::size_t >& addedTets,
const std::set< std::size_t >& removedNodes,
const std::unordered_map< int, std::unordered_set< std::size_t > >&
nodeCommMap,
const std::map< int, std::vector< std::size_t > >& bface,
const std::map< int, std::vector< std::size_t > >& bnode,
const std::vector< std::size_t >& triinpoel );
//! Const-ref access to current solution
//! \return Const-ref to current solution
const tk::Fields& solution() const { return m_u; }
//! Compute integral quantities for output
void integrals();
//! Compute recent conjugate gradients solution gradient
void sgrad();
//! Evaluate whether to continue with next time step
void step();
// Evaluate whether to do load balancing
void evalLB( int nrestart );
/** @name Charm++ pack/unpack serializer member functions */
///@{
//! \brief Pack/Unpack serialize member function
//! \param[in,out] p Charm++'s PUP::er serializer object reference
void pup( PUP::er &p ) override {
p | m_disc;
p | m_cgpre;
p | m_cgmom;
p | m_nrhs;
p | m_nnorm;
p | m_nsgrad;
p | m_npgrad;
p | m_nvgrad;
p | m_nflux;
p | m_ndiv;
p | m_nbpint;
p | m_np;
p | m_bnode;
p | m_bface;
p | m_triinpoel;
p | m_u;
p | m_un;
p | m_pr;
p | m_rhs;
p | m_rhsc;
p | m_sgrad;
p | m_sgradc;
p | m_pgrad;
p | m_pgradc;
p | m_vgrad;
p | m_vgradc;
p | m_flux;
p | m_fluxc;
p | m_div;
p | m_divc;
p | m_diag;
p | m_bnorm;
p | m_bnormc;
p | m_bndpoinint;
p | m_domedgeint;
p | m_bpint;
p | m_dsupedge;
p | m_dsupint;
p | m_dirbcmask;
p | m_dirbcval;
p | m_dirbcmaskp;
p | m_dirbcvalp;
p | m_symbcnodes;
p | m_symbcnorms;
p | m_noslipbcnodes;
p | m_surfint;
p | m_stage;
p | m_finished;
p | m_freezeflow;
p | m_rkcoef;
p | m_timer;
}
//! \brief Pack/Unpack serialize operator|
//! \param[in,out] p Charm++'s PUP::er serializer object reference
//! \param[in,out] i ChoCG object reference
friend void operator|( PUP::er& p, ChoCG& i ) { i.pup(p); }
//@}
private:
//! Discretization proxy
CProxy_Discretization m_disc;
//! Conjugate Gradients Charm++ proxy for pressure solve
tk::CProxy_ConjugateGradients m_cgpre;
//! Conjugate Gradients Charm++ proxy for momentum solve
tk::CProxy_ConjugateGradients m_cgmom;
//! Counter for right-hand side vector nodes updated
std::size_t m_nrhs;
//! Counter for receiving boundary point normals
std::size_t m_nnorm;
//! Counter for receiving conjugrate gradient solution gradient
std::size_t m_nsgrad;
//! Counter for receiving pressure gradient
std::size_t m_npgrad;
//! Counter for receiving velocity gradient
std::size_t m_nvgrad;
//! Counter for receiving momentum flux
std::size_t m_nflux;
//! Counter for receiving boundary velocity divergences
std::size_t m_ndiv;
//! Counter for receiving boundary point integrals
std::size_t m_nbpint;
//! Count number of Poisson solves during setup
std::size_t m_np;
//! Boundary node lists mapped to side set ids used in the input file
std::map< int, std::vector< std::size_t > > m_bnode;
//! Boundary face lists mapped to side set ids used in the input file
std::map< int, std::vector< std::size_t > > m_bface;
//! Boundary triangle face connecitivity where BCs are set by user
std::vector< std::size_t > m_triinpoel;
//! Unknown/solution vector at mesh nodes
tk::Fields m_u;
//! Unknown/solution vector at mesh nodes at previous time step
tk::Fields m_un;
//! Pressure
std::vector< tk::real > m_pr;
//! Right-hand side vector (for the high order system)
tk::Fields m_rhs;
//! Receive buffer for communication of the right hand side
//! \details Key: global node id, value: rhs for all scalar components per
//! node.
std::unordered_map< std::size_t, std::vector< tk::real > > m_rhsc;
//! Conjugate gradient solution gradient in mesh nodes
tk::Fields m_sgrad;
//! Conjugate gradient solution gradient receive buffer
std::unordered_map< std::size_t, std::vector< tk::real > > m_sgradc;
//! Pressure gradient in mesh nodes
tk::Fields m_pgrad;
//! Pressure gradient receive buffer
std::unordered_map< std::size_t, std::vector< tk::real > > m_pgradc;
//! Velocity gradient in mesh nodes
tk::Fields m_vgrad;
//! Velocity gradient receive buffer
std::unordered_map< std::size_t, std::vector< tk::real > > m_vgradc;
//! Momentum flux in mesh nodes
tk::Fields m_flux;
//! Momentum flux receive buffer
std::unordered_map< std::size_t, std::vector< tk::real > > m_fluxc;
//! Velocity divergence
std::vector< tk::real > m_div;
//! Receive buffer for communication of the velocity divergence
//! \details Key: global node id, value: velocity divergence
std::unordered_map< std::size_t, tk::real > m_divc;
//! Diagnostics object
NodeDiagnostics m_diag;
//! Boundary point normals
//! \details Outer key: side set id. Inner key: global node id of boundary
//! point, value: weighted normal vector, inverse distance square.
std::unordered_map< int,
std::unordered_map< std::size_t, std::array< tk::real, 4 > > > m_bnorm;
//! Boundary point normals receive buffer
//! \details Outer key: side set id. Inner key: global node id of boundary
//! point, value: weighted normals and inverse distance square.
decltype(m_bnorm) m_bnormc;
//! Boundary point integrals
//! \details Key: global node id of boundary point, value: boundary point
//! integral contributions.
std::unordered_map< std::size_t, std::array<tk::real,3> > m_bndpoinint;
//! Domain edge integrals
std::unordered_map< tk::UnsMesh::Edge, std::array< tk::real, 5 >,
tk::UnsMesh::Hash<2>, tk::UnsMesh::Eq<2> > m_domedgeint;
//! Streamable boundary point integrals
std::vector< tk::real > m_bpint;
//! Superedge (tet, face, edge) end points with local ids for domain edges
std::array< std::vector< std::size_t >, 3 > m_dsupedge;
//! Superedge (tet, face, edge) domain edge integrals
std::array< std::vector< tk::real >, 3 > m_dsupint;
//! Nodes and their Dirichlet BC masks
std::vector< std::size_t > m_dirbcmask;
//! Nodes and their Dirichlet BC values
std::vector< double > m_dirbcval;
//! Nodes and their pressure Dirichlet BC masks
std::vector< std::size_t > m_dirbcmaskp;
//! Nodes and their pressure Dirichlet BC values
std::vector< double > m_dirbcvalp;
//! Streamable nodes at which symmetry BCs are set
std::vector< std::size_t > m_symbcnodes;
//! Streamable normals at nodes at which symmetry BCs are set
std::vector< tk::real > m_symbcnorms;
//! Streamable nodes at which noslip BCs are set
std::vector< std::size_t > m_noslipbcnodes;
//! Streamable surface integral nodes and normals * dA on surfaces
std::map< int, std::pair< std::vector< std::size_t >,
std::vector< tk::real > > > m_surfint;
//! Runge-Kutta stage counter
std::size_t m_stage;
//! True in the last time step
int m_finished;
//! dt multiplier after flow no longer updated
tk::real m_freezeflow;
//! Runge-Kutta coefficients
std::vector< tk::real > m_rkcoef;
//! Timer
std::vector< tk::Timer > m_timer;
//! Access bound Discretization class pointer
Discretization* Disc() const {
Assert( m_disc[ thisIndex ].ckLocal() != nullptr, "ckLocal() null" );
return m_disc[ thisIndex ].ckLocal();
}
//! Access bound momentum solve matrix
tk::CSR& Lhs() {
Assert( m_cgmom[ thisIndex ].ckLocal() != nullptr, "ckLocal() null" );
return m_cgmom[ thisIndex ].ckLocal()->lhs();
}
//! Prepare Dirichlet boundary condition data structures
void setupDirBC( const std::vector< std::vector< int > >& cfgmask,
const std::vector< std::vector< double > >& cfgval,
std::size_t ncomp,
std::vector< std::size_t >& mask,
std::vector< double >& val );
//! Start computing velocity divergence
void div( const tk::Fields& u );
//! Start computing velocity gradient
void velgrad();
//! Start computing momentum flux
void flux();
//! Finalize computing gradient
void fingrad( tk::Fields& grad,
std::unordered_map< std::size_t, std::vector< tk::real > >& gradc );
//! Finalize computing pressure gradient
void finpgrad();
//! Compute pressure gradient
void pgrad();
//! Compute local contributions to domain edge integrals
void domint();
//! Setup lhs matrix for pressure solve
std::tuple< tk::CSR, std::vector< tk::real >, std::vector< tk::real > >
prelhs( const std::pair< std::vector< std::size_t >,
std::vector< std::size_t > >& psup );
//! Setup empty lhs matrix for momentum solve
std::tuple< tk::CSR, std::vector< tk::real >, std::vector< tk::real > >
momlhs( const std::pair< std::vector< std::size_t >,
std::vector< std::size_t > >& psup );
//! Compute chare-boundary edges
void bndEdges();
//! Compute local contributions to boundary normals and integrals
void bndint();
//! Combine own and communicated portions of the boundary point normals
void bnorm();
//! Convert integrals into streamable data structures
void streamable();
//! Generate superedge-groups for domain-edge loops
void domsuped();
//! Output mesh and particle fields to files
void out();
//! Output mesh-based fields to file
void writeFields( CkCallback cb );
//! Combine own and communicated portions of the integrals
void merge();
//! Fill lhs matrix of transport equations
void lhs();
//! Advance systems of equations
void solve();
//! Compute advective-diffusive prediction of momentum/transport
void pred();
//! Compute pressure correction
void corr();
//! Optionally refine/derefine mesh
void refine();
//! Compute time step size
void dt();
//! Evaluate whether to save checkpoint/restart
void evalRestart();
//! Apply boundary conditions
void BC( tk::Fields& u, tk::real t );
//! Apply scalar source to solution
void src();
};
} // inciter::
|