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

Change loops to use int instead of long #438

Merged
merged 1 commit into from
Apr 11, 2024
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
18 changes: 9 additions & 9 deletions Source/MaestroMakeEdgeState.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -60,7 +60,7 @@ void Maestro::MakeEdgeState1dSphr(BaseState<Real>& s_state,
const auto w0_arr = w0.const_array();

if (ppm_type == 0) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(int r) {
Real slope = 0.0;

// this will hold values at r-1, r and r+1
Expand Down Expand Up @@ -119,7 +119,7 @@ void Maestro::MakeEdgeState1dSphr(BaseState<Real>& s_state,
Gpu::synchronize();

} else if (ppm_type == 1) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(int r) {
// interpolate s to radial edges

// sm
Expand Down Expand Up @@ -220,7 +220,7 @@ void Maestro::MakeEdgeState1dSphr(BaseState<Real>& s_state,
});
Gpu::synchronize();
} else if (ppm_type == 2) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_fine, [=] AMREX_GPU_DEVICE(int r) {
// interpolate s to radial edges, store these temporary values into sedgel

// this will hold values at r-1, r, r+1 and r+2
Expand Down Expand Up @@ -364,7 +364,7 @@ void Maestro::MakeEdgeState1dSphr(BaseState<Real>& s_state,
Gpu::synchronize();
}

ParallelFor(nr_fine + 1, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_fine + 1, [=] AMREX_GPU_DEVICE(int r) {
// Fix center and edge of star by reflecting the extrapolated state.
// An alternate way would be to compute these values using the entire algorithm,
// but that would require more ghost cells at several stages.
Expand Down Expand Up @@ -439,7 +439,7 @@ void Maestro::MakeEdgeState1dPlanar(BaseState<Real>& s_state,

if (ppm_type == 0) {
// compute slopes
ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(long j) {
ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(int j) {
Real slope = 0.0;
int r = j + lo;

Expand Down Expand Up @@ -563,7 +563,7 @@ void Maestro::MakeEdgeState1dPlanar(BaseState<Real>& s_state,
} else if (ppm_type == 1) {
// interpolate s to radial edges, store these temporary values into sedgel

ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(long j) {
ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(int j) {
int r = j + lo;

// calculate sm
Expand Down Expand Up @@ -718,7 +718,7 @@ void Maestro::MakeEdgeState1dPlanar(BaseState<Real>& s_state,
BaseState<Real> sedget_s(base_geom.nr_fine + 1);
auto sedget = sedget_s.array();

ParallelFor(hi - lo + 2, [=] AMREX_GPU_DEVICE(long j) {
ParallelFor(hi - lo + 2, [=] AMREX_GPU_DEVICE(int j) {
int r = j + lo;

// left side
Expand Down Expand Up @@ -784,7 +784,7 @@ void Maestro::MakeEdgeState1dPlanar(BaseState<Real>& s_state,
});
Gpu::synchronize();

ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(long j) {
ParallelFor(hi - lo + 1, [=] AMREX_GPU_DEVICE(int j) {
int r = j + lo;

// use Colella 2008 limiters
Expand Down Expand Up @@ -957,7 +957,7 @@ void Maestro::MakeEdgeState1dPlanar(BaseState<Real>& s_state,
const int hi = base_geom.r_end_coord(n, i);

// solve Riemann problem to get final edge state
ParallelFor(hi - lo + 2, [=] AMREX_GPU_DEVICE(long j) {
ParallelFor(hi - lo + 2, [=] AMREX_GPU_DEVICE(int j) {
int r = j + lo;

if (r == 0) {
Expand Down
4 changes: 2 additions & 2 deletions Source/MaestroMakeGrav.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -21,7 +21,7 @@ void Maestro::MakeGravCell(BaseState<Real>& grav_cell,
// does not contribute to the gravitational acceleration.
for (auto n = 0; n <= base_geom.finest_radial_level; ++n) {
const int nr_lev = base_geom.nr(n);
ParallelFor(nr_lev, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_lev, [=] AMREX_GPU_DEVICE(int r) {
grav_cell_arr(n, r) = -Gconst * planar_invsq_mass_loc /
(r_cc_loc(n, r) * r_cc_loc(n, r));
});
Expand Down Expand Up @@ -223,7 +223,7 @@ void Maestro::MakeGravEdge(BaseState<Real>& grav_edge_state,
//
for (auto n = 0; n <= base_geom.finest_radial_level; ++n) {
const int nr_lev = base_geom.nr(n);
ParallelFor(nr_lev, [=] AMREX_GPU_DEVICE(long r) {
ParallelFor(nr_lev, [=] AMREX_GPU_DEVICE(int r) {
grav_edge(n, r) = -Gconst * planar_invsq_mass_loc /
(r_edge_loc(n, r) * r_edge_loc(n, r));
});
Expand Down