Skip to content
Snippets Groups Projects

Draft: C++ port of horizontal/mo_lib_gradients.F90

Open Gwenolé Lucas requested to merge feature-add-cpp-codes-gradients into feature-add-cpp-codes
2 files
+ 70
1
Compare changes
  • Side-by-side
  • Inline
Files
2
@@ -66,3 +66,60 @@ void grad_fd_norm_lib(const T* psi_c,
}
}
// Computes directional derivative of a vertex centered variable with
// respect to direction tangent to triangle edge. Notice that the
// tangential direction is defined by: iorient*(vertex2 - vertex1)
// input: lives on vertices of triangles
// output: lives on edges (velocity points)
template <typename T>
void grad_fd_norm_lib(const T* psi_v,
const int* edge_vertex_idx, const int* edge_vertex_blk,
const T* primal_edge_length, const T* tangent_orientation,
T* grad_tang_psi_e,
// Integer parameters
int i_startblk, int i_endblk,
int i_startidx_in, int i_endidx_in,
int slev, int elev, int nproma,
// Extra array dimension
int nlev, int nblks_v, int nblks_e,
// OpenACC flag
bool lacc)
{
// Host or device check
// Wrap raw pointers in unmanaged Kokkos views
typedef Kokkos::View<const T***, Kokkos::LayoutLeft, Kokkos::MemoryUnmanaged> UnmanagedConstT3D;
typedef Kokkos::View<const int***, Kokkos::LayoutLeft, Kokkos::MemoryUnmanaged> UnmanagedConstInt3D;
typedef Kokkos::View<const int**, Kokkos::LayoutLeft, Kokkos::MemoryUnmanaged> UnmanagedConstT2D;
typedef Kokkos::View<T***, Kokkos::LayoutLeft, Kokkos::MemoryUnmanaged> UnmanagedT3D;
UnmanagedConstT3D psi_v_view(psi_v, nproma, nlev, nblks_v);
UnmanagedConstInt3D edge_vertex_idx_view(edge_vertex_idx, nproma, nblks_e, 2);
UnmanagedConstInt3D edge_vertex_blk_view(edge_vertex_blk, nproma, nblks_e, 2);
UnmanagedConstT2D primal_edge_length_view(primal_edge_length, nproma, nblks_e);
UnmanagedConstT2D tangent_orientation_view(tangent_orientation, nproma, nblks_e);
UnmanagedT3D grad_tang_psi_e_view(grad_tang_psi_e, nproma, nlev, nblks_e);
for (int jb = i_startblk; jb < i_endblk; ++jb)
{
int i_startidx, i_endidx;
get_indices_e_lib(i_startidx_in, i_endidx_in, nproma, jb,
i_startblk, i_endblk, i_startidx, i_endidx);
Kokkos::MDRangePolicy<Kokkos::Rank<2>> innerPolicy({slev, i_startidx},
{elev, i_endidx});
Kokkos::parallel_for("grad_fd_norm_lib_inner", innerPolicy,
KOKKOS_LAMBDA(const int jk, const int je)
{
grad_tang_psi_e_view(je, jk, jb) =
tangent_orientation_view(je, jb) *
(psi_v_view(edge_vertex_idx_view(je, jb, 1), jk, edge_vertex_blk_view(je, jb, 1)) -
psi_v_view(edge_vertex_idx_view(je, jb, 0), jk, edge_vertex_blk_view(je, jb, 0))) /
primal_edge_length_view(je, jb);
});
// Is a barrier needed here?
// Kokkos::fence();
}
}
Loading