Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

NodalProjector: fix comments #1396

Merged
merged 2 commits into from
Sep 23, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
27 changes: 23 additions & 4 deletions Src/LinearSolvers/Projections/AMReX_NodalProjector.H
Original file line number Diff line number Diff line change
Expand Up @@ -7,22 +7,41 @@
#include <AMReX_MLNodeLaplacian.H>
#include <AMReX_MLMG.H>

//
//
// *************************** DEFAULT MODE ***************************
//
// Solves
//
// div(sigma*grad(phi)) = div(vel) + S_nd + S_cc
//
// and performs the projection
//
// vel = vel - ( sigma / alpha ) * grad(phi)
// vel = vel - sigma * grad(phi)
//
// vel, sigma, alpha, and S_cc are cell-centered variables, while
// vel, sigma, and S_cc are cell-centered variables, while
// phi and S_nd are nodal-centered variables.
//
// *************************** CUSTOM MODE ***************************
//
// Solves
//
// div(sigma*grad(phi)) = rhs
//
// and performs the projection
//
// vel = vel - (sigma/alpha) * grad(phi)
//
// alpha is a cell-centered variable, while rhs is nodal-centered.
//
// In this mode, the user provides rhs and alpha
//
// By default alpha is assumed to be 1. Use setAlpha to change the default.
//
// The user can provide a custom RHS instead of letting NodalProjector to compute
// one. Use setCustomRHS to provide a custom RHS.
// Use setCustomRHS to provide a custom RHS, else "div(vel) + S_nd + S_cc"
// is used.
//
// Example: rhs = div(alpha*vel)
//
namespace amrex {

Expand Down
6 changes: 3 additions & 3 deletions Src/LinearSolvers/Projections/AMReX_NodalProjector.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -241,7 +241,7 @@ NodalProjector::project ( Real a_rtol, Real a_atol )
// phi comes out already averaged-down and ready to be used by caller if needed
m_mlmg -> solve( GetVecOfPtrs(m_phi), GetVecOfConstPtrs(m_rhs), a_rtol, a_atol );

// Get fluxes -- fluxes = - (alpha/beta) * grad(phi)
// Get fluxes -- fluxes = - sigma * grad(phi)
m_mlmg -> getFluxes( GetVecOfPtrs(m_fluxes) );

// At this time, the fluxes are "correct" only on regions not covered by finer grids.
Expand All @@ -258,15 +258,15 @@ NodalProjector::project ( Real a_rtol, Real a_atol )
{
if (m_has_alpha)
{
// fluxes -> fluxes/alpha = -grad(phi)/beta
// fluxes -> fluxes/alpha = - ( sigma / alpha ) * grad(phi)
for (int n = 0; n < AMREX_SPACEDIM; ++n)
{
MultiFab::Divide( m_fluxes[lev], *m_alpha[lev], 0, n, 1, 0 );
}
}

//
// vel = vel + fluxes = vel - grad(phi) / beta
// vel = vel + fluxes = vel - ( sigma / alpha ) * grad(phi)
//
// Since we already averaged-down the velocity field and -grad(phi),
// we perform the projection by simply adding the two of them.
Expand Down