inciter::ChoCG class

ChoCG Charm++ chare array used to advance PDEs in time with ChoCG.

Public static functions

static void registerReducers()
Configure Charm++ custom reduction types initiated from this chare array.

Constructors, destructors, conversion operators

ChoCG(CkMigrateMessage* m) explicit
Migrate constructor.

Public functions

auto 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) -> ChoCG_SDAG_CODE explicit
Constructor.
void ResumeFromSync() override
Return from migration.
void setup(tk::real v)
Start setup for solution.
void pinit()
Initialize Poisson solve.
void psolve()
Solve Poisson equation.
void psolved()
Continue after Poisson solve.
void msolve()
Solve momentum/transport equations.
void msolved()
Continue after momentum/transport ssolve.
void advance(tk::real newdt)
Advance equations to next time step.
void diag()
Evaluate diagnostics.
void feop()
Start (re-)computing domain and boundary integrals.
void comnorm(const std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>>& inbnd)
Receive contributions to boundary point normals on chare-boundaries.
void comvgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)
Receive contributions to velocity gradients.
void comflux(const std::unordered_map<std::size_t, std::vector<tk::real>>& influx)
Receive contributions to momentum flux on chare-boundaries.
void comsgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)
Receive contributions to conjugate gradients solution gradient.
void compgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)
Receive contributions to pressure gradient.
void comrhs(const std::unordered_map<std::size_t, std::vector<tk::real>>& inrhs)
Receive contributions to right-hand side vector on chare-boundaries.
void comdiv(const std::unordered_map<std::size_t, tk::real>& indiv)
Receive contributions to velocity divergence on chare-boundaries.
void comaec(const std::unordered_map<std::size_t, std::vector<tk::real>>& inaec)
Receive antidiffusive and low-order contributions on chare-boundaries.
void comalw(const std::unordered_map<std::size_t, std::vector<tk::real>>& inalw)
Receive allowed limits contributions on chare-boundaries.
void comlim(const std::unordered_map<std::size_t, std::vector<tk::real>>& inlim)
Receive limited antidiffusive contributions on chare-boundaries.
void evalres(const std::vector<tk::real>& l2res)
Evaluate residuals.
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)
Receive new mesh from Refiner.
auto solution() const -> const tk::Fields&
void integrals()
Compute integral quantities for output.
void sgrad()
Compute recent conjugate gradients solution gradient.
void step()
Evaluate whether to continue with next time step.
void evalLB(int nrestart)

Charm++ pack/unpack serializer member functions

CProxy_Discretization m_disc
Discretization proxy.
tk::CProxy_ConjugateGradients m_cgpre
Conjugate Gradients Charm++ proxy for pressure solve.
tk::CProxy_ConjugateGradients m_cgmom
Conjugate Gradients Charm++ proxy for momentum solve.
std::size_t m_nrhs
Counter for right-hand side vector nodes updated.
std::size_t m_nnorm
Counter for receiving boundary point normals.
std::size_t m_naec
Counter for receiving antidiffusive contributions.
std::size_t m_nalw
Counter for receiving allowed limits.
std::size_t m_nlim
Counter for receiving limited antidiffusive contributions.
std::size_t m_nsgrad
Counter for receiving conjugrate gradient solution gradient.
std::size_t m_npgrad
Counter for receiving pressure gradient.
std::size_t m_nvgrad
Counter for receiving velocity gradient.
std::size_t m_nflux
Counter for receiving momentum flux.
std::size_t m_ndiv
Counter for receiving boundary velocity divergences.
std::size_t m_nbpint
Counter for receiving boundary point integrals.
std::size_t m_np
Count number of Poisson solves during setup.
std::map<int, std::vector<std::size_t>> m_bnode
Boundary node lists mapped to side set ids used in the input file.
std::map<int, std::vector<std::size_t>> m_bface
Boundary face lists mapped to side set ids used in the input file.
std::vector<std::size_t> m_triinpoel
Boundary triangle face connecitivity where BCs are set by user.
tk::Fields m_u
Unknown/solution vector at mesh nodes.
tk::Fields m_un
Unknown/solution vector at mesh nodes at previous time step.
std::vector<tk::real> m_pr
Pressure.
tk::Fields m_p
Max/min antidiffusive edge contributions at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_pc
tk::Fields m_q
Max/min allowed limits at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_qc
tk::Fields m_a
Limited antidiffusive contributions at mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_ac
tk::Fields m_rhs
Right-hand side vector (for the high order system)
std::unordered_map<std::size_t, std::vector<tk::real>> m_rhsc
tk::Fields m_sgrad
Conjugate gradient solution gradient in mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_sgradc
Conjugate gradient solution gradient receive buffer.
tk::Fields m_pgrad
Pressure gradient in mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_pgradc
Pressure gradient receive buffer.
tk::Fields m_vgrad
Velocity gradient in mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_vgradc
Velocity gradient receive buffer.
tk::Fields m_flux
Momentum flux in mesh nodes.
std::unordered_map<std::size_t, std::vector<tk::real>> m_fluxc
Momentum flux receive buffer.
std::vector<tk::real> m_div
Velocity divergence.
std::unordered_map<std::size_t, tk::real> m_divc
NodeDiagnostics m_diag
Diagnostics object.
std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>> m_bnorm
decltype(m_bnorm) m_bnormc
std::unordered_map<std::size_t, std::array<tk::real, 3>> m_bndpoinint
std::unordered_map<tk::UnsMesh::Edge, std::array<tk::real, 5>, tk::UnsMesh::Hash<2>, tk::UnsMesh::Eq<2>> m_domedgeint
Domain edge integrals.
std::vector<tk::real> m_bpint
Streamable boundary point integrals.
std::array<std::vector<std::size_t>, 3> m_dsupedge
Superedge (tet, face, edge) end points with local ids for domain edges.
std::array<std::vector<tk::real>, 3> m_dsupint
Superedge (tet, face, edge) domain edge integrals.
std::array<std::vector<tk::real>, 3> m_dsuplim
FCT limiter coefficients in domain superedges.
std::vector<std::size_t> m_dirbcmask
Nodes and their Dirichlet BC masks.
std::vector<double> m_dirbcval
Nodes and their Dirichlet BC values.
std::vector<std::size_t> m_dirbcmaskp
Nodes and their pressure Dirichlet BC masks.
std::vector<double> m_dirbcvalp
Nodes and their pressure Dirichlet BC values.
std::vector<std::size_t> m_symbcnodes
Streamable nodes at which symmetry BCs are set.
std::vector<tk::real> m_symbcnorms
Streamable normals at nodes at which symmetry BCs are set.
std::vector<std::size_t> m_noslipbcnodes
Streamable nodes at which noslip BCs are set.
std::map<int, std::pair<std::vector<std::size_t>, std::vector<tk::real>>> m_surfint
Streamable surface integral nodes and normals * dA on surfaces.
std::size_t m_stage
Runge-Kutta stage counter.
int m_finished
True in the last time step.
std::vector<tk::real> m_rkcoef
Runge-Kutta coefficients.
std::vector<tk::Timer> m_timer
Timer.
void pup(PUP::er& p) override
Pack/Unpack serialize member function.
void operator|(PUP::er& p, ChoCG& i)
Pack/Unpack serialize operator|.

Function documentation

static void inciter::ChoCG::registerReducers()

Configure Charm++ custom reduction types initiated from this chare array.

Since this is a [initnode] routine, the runtime system executes the routine exactly once on every logical node early on in the Charm++ init sequence. Must be static as it is called without an object. See also: Section "Initializations at Program Startup" at in the Charm++ manual http://charm.cs.illinois.edu/manuals/html/charm++/manual.html.

ChoCG_SDAG_CODE inciter::ChoCG::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) explicit

Constructor.

Parameters
disc in Discretization proxy
cgpre in ConjugateGradients Charm++ proxy for pressure solve
cgmom in ConjugateGradients Charm++ proxy for momentum solve
bface in Boundary-faces mapped to side sets used in the input file
bnode in Boundary-node lists mapped to side sets used in input file
triinpoel in Boundary-face connectivity where BCs set (global ids)

void inciter::ChoCG::ResumeFromSync() override

Return from migration.

This is called when load balancing (LB) completes. The presence of this function does not affect whether or not we block on LB.

void inciter::ChoCG::setup(tk::real v)

Start setup for solution.

Parameters
in Total volume within user-specified box

void inciter::ChoCG::advance(tk::real newdt)

Advance equations to next time step.

Parameters
newdt in The smallest dt across the whole problem

void inciter::ChoCG::comnorm(const std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>>& inbnd)

Receive contributions to boundary point normals on chare-boundaries.

Parameters
inbnd in Incoming partial sums of boundary point normals

void inciter::ChoCG::comvgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)

Receive contributions to velocity gradients.

Parameters
ingrad in Partial contributions of momentum flux to chare-boundary nodes. Key: global mesh node IDs, values: contributions.

This function receives contributions to m_vgrad, which stores the velocity gradients at mesh nodes. While m_vgrad stores own contributions, m_vgradc collects the neighbor chare contributions during communication. This way work on m_vgrad and m_vgradc is overlapped.

void inciter::ChoCG::comflux(const std::unordered_map<std::size_t, std::vector<tk::real>>& influx)

Receive contributions to momentum flux on chare-boundaries.

Parameters
influx in Partial contributions of momentum flux to chare-boundary nodes. Key: global mesh node IDs, values: contributions.

This function receives contributions to m_flux, which stores the momentum flux at mesh nodes. While m_flux stores own contributions, m_fluxc collects the neighbor chare contributions during communication. This way work on m_flux and m_fluxc is overlapped.

void inciter::ChoCG::comsgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)

Receive contributions to conjugate gradients solution gradient.

Parameters
ingrad in Partial contributions to chare-boundary nodes. Key: global mesh node IDs, value: contributions for all scalar components.

This function receives contributions to m_sgrad, which stores the gradients at mesh nodes. While m_sgrad stores own contributions, m_sgradc collects the neighbor chare contributions during communication. This way work on m_sgrad and m_sgradc is overlapped.

void inciter::ChoCG::compgrad(const std::unordered_map<std::size_t, std::vector<tk::real>>& ingrad)

Receive contributions to pressure gradient.

Parameters
ingrad in Partial contributions to chare-boundary nodes. Key: global mesh node IDs, value: contributions for all scalar components.

This function receives contributions to m_pgrad, which stores the gradients at mesh nodes. While m_pgrad stores own contributions, m_pgradc collects the neighbor chare contributions during communication. This way work on m_pgrad and m_pgradc is overlapped.

void inciter::ChoCG::comrhs(const std::unordered_map<std::size_t, std::vector<tk::real>>& inrhs)

Receive contributions to right-hand side vector on chare-boundaries.

Parameters
inrhs in Partial contributions of RHS to chare-boundary nodes. Key: global mesh node IDs, value: contributions for all scalar components.

This function receives contributions to m_rhs, which stores the right hand side vector at mesh nodes. While m_rhs stores own contributions, m_rhsc collects the neighbor chare contributions during communication. This way work on m_rhs and m_rhsc is overlapped.

void inciter::ChoCG::comdiv(const std::unordered_map<std::size_t, tk::real>& indiv)

Receive contributions to velocity divergence on chare-boundaries.

Parameters
indiv in Partial contributions of velocity divergence to chare-boundary nodes. Key: global mesh node IDs, value: contribution.

This function receives contributions to m_div, which stores the velocity divergence at mesh nodes. While m_div stores own contributions, m_divc collects the neighbor chare contributions during communication. This way work on m_div and m_divc is overlapped.

void inciter::ChoCG::comaec(const std::unordered_map<std::size_t, std::vector<tk::real>>& inaec)

Receive antidiffusive and low-order contributions on chare-boundaries.

Parameters
inaec in Partial contributions of antidiffusive edge and low-order solution contributions on chare-boundary nodes. Key: global mesh node IDs, value: 0: antidiffusive contributions, 1: low-order solution.

void inciter::ChoCG::comalw(const std::unordered_map<std::size_t, std::vector<tk::real>>& inalw)

Receive allowed limits contributions on chare-boundaries.

Parameters
inalw in Partial contributions of allowed limits contributions on chare-boundary nodes. Key: global mesh node IDs, value: allowed limit contributions.

void inciter::ChoCG::comlim(const std::unordered_map<std::size_t, std::vector<tk::real>>& inlim)

Receive limited antidiffusive contributions on chare-boundaries.

Parameters
inlim in Partial contributions of limited contributions on chare-boundary nodes. Key: global mesh node IDs, value: limited contributions.

void inciter::ChoCG::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)

Receive new mesh from Refiner.

Parameters
ginpoel in Mesh connectivity with global node ids
chunk in New mesh chunk (connectivity and global<->local id maps)
coord in New mesh node coordinates
addedNodes in Newly added mesh nodes and their parents (local ids)
addedTets in Newly added mesh cells and their parents (local ids)
removedNodes in Newly removed mesh node local ids
nodeCommMap in New node communication map
bface in Boundary-faces mapped to side set ids
bnode in Boundary-node lists mapped to side set ids
triinpoel in Boundary-face connectivity

const tk::Fields& inciter::ChoCG::solution() const

Returns Const-ref to current solution

Const-ref access to current solution

void inciter::ChoCG::evalLB(int nrestart)

Parameters
nrestart in Number of times restarted

void inciter::ChoCG::pup(PUP::er& p) override

Pack/Unpack serialize member function.

Parameters
in/out Charm++'s PUP::er serializer object reference

void inciter::ChoCG::operator|(PUP::er& p, ChoCG& i)

Pack/Unpack serialize operator|.

Parameters
in/out Charm++'s PUP::er serializer object reference
in/out ChoCG object reference

Variable documentation

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::ChoCG::m_pc

Receive buffer for max/min antidiffusive edge contributions

Key: global node id, value: max/min antidiff edge contributions

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::ChoCG::m_qc

Receive buffer for max/min allowed limits

Key: global node id, value: max/min allowed limits in nodes.

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::ChoCG::m_ac

Receive buffer for limited antidiffusive contributions

Key: global node id, value: limited antidiffusive contributions in nodes.

std::unordered_map<std::size_t, std::vector<tk::real>> inciter::ChoCG::m_rhsc

Receive buffer for communication of the right hand side

Key: global node id, value: rhs for all scalar components per node.

std::unordered_map<std::size_t, tk::real> inciter::ChoCG::m_divc

Receive buffer for communication of the velocity divergence

Key: global node id, value: velocity divergence

std::unordered_map<int, std::unordered_map<std::size_t, std::array<tk::real, 4>>> inciter::ChoCG::m_bnorm

Boundary point normals

Outer key: side set id. Inner key: global node id of boundary point, value: weighted normal vector, inverse distance square.

decltype(m_bnorm) inciter::ChoCG::m_bnormc

Boundary point normals receive buffer

Outer key: side set id. Inner key: global node id of boundary point, value: weighted normals and inverse distance square.

std::unordered_map<std::size_t, std::array<tk::real, 3>> inciter::ChoCG::m_bndpoinint

Boundary point integrals

Key: global node id of boundary point, value: boundary point integral contributions.