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

Update profile output for non-zero center points #100

Open
wants to merge 2 commits into
base: main
Choose a base branch
from
Open
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
36 changes: 29 additions & 7 deletions src/CabanaPD_DisplacementProfile.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -42,30 +42,52 @@ void createOutputProfile( MPI_Comm comm, const int num_cell,
Kokkos::View<int*, memory_space> count( "c", 1 );

auto dims = getDim( profile_dim );
double dx1 = particles.dx[dims[0]];
double dx2 = particles.dx[dims[1]];
double dx1 = particles.dx[dims[0]] / 2.0;
double dx2 = particles.dx[dims[1]] / 2.0;
double center1 = particles.global_mesh_lo[dims[0]] +
particles.global_mesh_ext[dims[0]] / 2.0;
double center2 = particles.global_mesh_lo[dims[1]] +
particles.global_mesh_ext[dims[1]] / 2.0;

auto x = particles.sliceReferencePosition();
auto u = particles.sliceDisplacement();
// Find points closest to the center line.
double center_min1;
double center_min2;
auto find_profile =
KOKKOS_LAMBDA( const int pid, double& min1, double& min2 )
{
if ( Kokkos::abs( x( pid, dims[0] ) - center1 ) < min1 )
min1 = x( pid, dims[0] );
if ( Kokkos::abs( x( pid, dims[1] ) - center2 ) < min2 )
min2 = x( pid, dims[1] );
};
Kokkos::RangePolicy<typename memory_space::execution_space> policy(
0, x.size() );
Kokkos::parallel_reduce( "displacement_profile", policy, find_profile,
Kokkos::Min<double>( center_min1 ),
Kokkos::Min<double>( center_min2 ) );

// Extract points along that line.
auto measure_profile = KOKKOS_LAMBDA( const int pid )
{
if ( x( pid, dims[0] ) < dx1 / 2.0 && x( pid, dims[0] ) > -dx1 / 2.0 &&
x( pid, dims[1] ) < dx2 / 2.0 && x( pid, dims[1] ) > -dx2 / 2.0 )
if ( x( pid, dims[0] ) - center_min1 < dx1 &&
x( pid, dims[0] ) - center_min1 > -dx1 &&
x( pid, dims[1] ) - center_min2 < dx2 &&
x( pid, dims[1] ) - center_min2 > -dx2 )
{
auto c = Kokkos::atomic_fetch_add( &count( 0 ), 1 );
profile( c, 0 ) = x( pid, profile_dim );
profile( c, 1 ) = user( pid );
}
};
Kokkos::RangePolicy<typename memory_space::execution_space> policy(
0, x.size() );
Kokkos::parallel_for( "displacement_profile", policy, measure_profile );
auto count_host =
Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace{}, count );
auto profile_host =
Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace{}, profile );
std::fstream fout;

std::fstream fout;
fout.open( file_name, std::ios::app );
for ( int p = 0; p < count_host( 0 ); p++ )
{
Expand Down
6 changes: 6 additions & 0 deletions src/CabanaPD_Particles.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -125,6 +125,8 @@ class Particles<MemorySpace, PMB, Dimension>

// Simulation total domain.
std::array<double, dim> global_mesh_ext;
std::array<double, dim> global_mesh_lo;
std::array<double, dim> global_mesh_hi;

// Simulation sub domain (single MPI rank).
std::array<double, dim> local_mesh_ext;
Expand Down Expand Up @@ -201,6 +203,8 @@ class Particles<MemorySpace, PMB, Dimension>
std::array<bool, dim> is_periodic;
for ( int d = 0; d < dim; d++ )
{
global_mesh_lo[d] = global_mesh->lowCorner( d );
global_mesh_hi[d] = global_mesh->highCorner( d );
global_mesh_ext[d] = global_mesh->extent( d );
is_periodic[d] = false;
}
Expand Down Expand Up @@ -516,6 +520,8 @@ class Particles<MemorySpace, LPS, Dimension>

// Simulation total domain.
using base_type::global_mesh_ext;
using base_type::global_mesh_hi;
using base_type::global_mesh_lo;

// Simulation sub domain (single MPI rank).
using base_type::ghost_mesh_hi;
Expand Down
Loading