4#include "../../quadrature/quadrature.hpp"
16template <
typename ScalarT >
17KOKKOS_INLINE_FUNCTION ScalarT
18 supg_tau(
const ScalarT vel_norm,
const ScalarT kappa,
const ScalarT h,
const ScalarT Pe_tol = 1e-8 )
20 const ScalarT kappa_min = ScalarT( 1e-8 );
21 const ScalarT kappa_virtual = Kokkos::max( kappa, kappa_min );
24 const ScalarT eps_vel = ScalarT( 1e-12 );
26 if ( vel_norm <= eps_vel )
32 const ScalarT Pe = vel_norm * h / ( 2.0 * kappa_virtual );
40 const ScalarT tau_max = ScalarT( 10.0 ) * h / Kokkos::max( vel_norm, eps_vel );
44 const ScalarT coth_term = ScalarT( 1.0 ) / Kokkos::tanh( Pe ) - ScalarT( 1.0 ) / Pe;
45 const ScalarT tau = ( h / ( 2.0 * vel_norm ) ) * coth_term;
47 return ( tau > tau_max ) ? tau_max : tau;
129template <
typename ScalarT,
int VelocityVecDim = 3 >
146 ScalarT diffusivity_;
149 bool treat_boundary_;
151 ScalarT mass_scaling_;
171 const ScalarT diffusivity,
174 bool diagonal =
false,
175 ScalarT mass_scaling = 1.0,
176 bool lumped_mass =
false,
183 , boundary_mask_( boundary_mask )
184 , velocity_( velocity )
185 , diffusivity_( diffusivity )
187 , treat_boundary_( treat_boundary )
188 , diagonal_( diagonal )
189 , mass_scaling_( mass_scaling )
190 , lumped_mass_( lumped_mass )
191 , operator_apply_mode_( operator_apply_mode )
192 , operator_communication_mode_( operator_communication_mode )
194 , send_buffers_( domain )
195 , recv_buffers_( domain )
198 ScalarT&
dt() {
return dt_; }
199 const ScalarT&
dt()
const {
return dt_; }
224 domain_, dst_, send_buffers_, recv_buffers_ );
229 KOKKOS_INLINE_FUNCTION
void
230 operator()(
const int local_subdomain_id,
const int x_cell,
const int y_cell,
const int r_cell )
const
237 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
238 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
244 ScalarT quad_weights[num_quad_points];
254 for (
int d = 0; d < VelocityVecDim; d++ )
257 vel_coeffs[d], local_subdomain_id, x_cell, y_cell, r_cell, d, vel_grid_ );
262 for (
int q = 0; q < num_quad_points; q++ )
266 const auto shape_i =
shape( i, quad_points[q] );
267 for (
int d = 0; d < VelocityVecDim; d++ )
269 vel_interp[wedge][q]( d ) += vel_coeffs[d][wedge]( i ) * shape_i;
280 const auto h = r_2 - r_1;
284 ScalarT tau_accum = 0.0;
285 ScalarT waccum = 0.0;
287 for (
int q = 0; q < num_quad_points; ++q )
290 const auto& uq = vel_interp[wedge][q];
291 const ScalarT vel_norm_q = uq.
norm();
293 const ScalarT tau_q = supg_tau< ScalarT >( vel_norm_q, diffusivity_, h, 1e-08 );
296 const ScalarT wq = quad_weights[q];
297 tau_accum += tau_q * wq;
302 ScalarT tau_cell = ( waccum > 0.0 ) ? ( tau_accum / waccum ) : 0.0;
303 streamline_diffusivity[wedge] = tau_cell;
310 for (
int q = 0; q < num_quad_points; q++ )
312 const auto w = quad_weights[q];
316 const auto J =
jac( wedge_phy_surf[wedge], r_1, r_2, quad_points[q] );
317 const auto det = Kokkos::abs( J.det() );
318 const auto J_inv_transposed = J.inv().transposed();
320 const auto vel = vel_interp[wedge][q];
324 const auto shape_i =
shape( i, quad_points[q] );
325 const auto grad_i = J_inv_transposed *
grad_shape( i, quad_points[q] );
329 const auto shape_j =
shape( j, quad_points[q] );
330 const auto grad_j = J_inv_transposed *
grad_shape( j, quad_points[q] );
332 const auto mass = shape_i * shape_j;
333 const auto diffusion = diffusivity_ * ( grad_i ).dot( grad_j );
334 const auto advection = ( vel.dot( grad_j ) ) * shape_i;
335 const auto streamline_diffusion =
336 streamline_diffusivity[wedge] * ( vel.dot( grad_j ) ) * ( vel.dot( grad_i ) );
338 M[wedge]( i, j ) += w * mass_scaling_ * mass * det;
339 A[wedge]( i, j ) += w * dt_ * ( diffusion + advection + streamline_diffusion ) * det;
353 if ( treat_boundary_ )
358 boundary_mask_( local_subdomain_id, x_cell, y_cell, r_cell + 1 ),
364 boundary_mask.
fill( 1.0 );
365 if ( at_cmb_boundary )
368 for (
int i = 0; i < 6; i++ )
370 for (
int j = 0; j < 6; j++ )
372 if ( i != j && ( i < 3 || j < 3 ) )
374 boundary_mask( i, j ) = 0.0;
380 if ( at_surface_boundary )
383 for (
int i = 0; i < 6; i++ )
385 for (
int j = 0; j < 6; j++ )
387 if ( i != j && ( i >= 3 || j >= 3 ) )
389 boundary_mask( i, j ) = 0.0;
413 dst[0] = ( M[0] + A[0] ) * src[0];
414 dst[1] = ( M[1] + A[1] ) * src[1];
Linear operator for a method-of-lines discretization of the unsteady advection-diffusion equation wit...
Definition unsteady_advection_diffusion_supg.hpp:131
UnsteadyAdvectionDiffusionSUPG(const grid::shell::DistributedDomain &domain, const grid::Grid3DDataVec< ScalarT, 3 > &grid, const grid::Grid2DDataScalar< ScalarT > &radii, const grid::Grid4DDataScalar< grid::shell::ShellBoundaryFlag > &boundary_mask, const linalg::VectorQ1Vec< ScalarT, VelocityVecDim > &velocity, const ScalarT diffusivity, const ScalarT dt, bool treat_boundary, bool diagonal=false, ScalarT mass_scaling=1.0, bool lumped_mass=false, linalg::OperatorApplyMode operator_apply_mode=linalg::OperatorApplyMode::Replace, linalg::OperatorCommunicationMode operator_communication_mode=linalg::OperatorCommunicationMode::CommunicateAdditively)
Definition unsteady_advection_diffusion_supg.hpp:165
void apply_impl(const SrcVectorType &src, DstVectorType &dst)
Definition unsteady_advection_diffusion_supg.hpp:201
void operator()(const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell) const
Definition unsteady_advection_diffusion_supg.hpp:230
const ScalarT & dt() const
Definition unsteady_advection_diffusion_supg.hpp:199
ScalarT ScalarType
Definition unsteady_advection_diffusion_supg.hpp:135
ScalarT & dt()
Definition unsteady_advection_diffusion_supg.hpp:198
Parallel data structure organizing the thick spherical shell metadata for distributed (MPI parallel) ...
Definition spherical_shell.hpp:2498
Q1 scalar finite element vector on a distributed shell grid.
Definition vector_q1.hpp:21
const grid::Grid4DDataScalar< ScalarType > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:137
Static assertion: VectorQ1Scalar satisfies VectorLike concept.
Definition vector_q1.hpp:162
const grid::Grid4DDataVec< ScalarType, VecDim > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:280
Timer supporting RAII scope or manual stop.
Definition timer.hpp:270
void stop()
Stop the timer and record elapsed time.
Definition timer.hpp:289
Concept for types that behave like linear operators.
Definition operator.hpp:57
void unpack_and_reduce_local_subdomain_boundaries(const grid::shell::DistributedDomain &domain, const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_recv_buffers, CommunicationReduction reduction=CommunicationReduction::SUM)
Unpacks and reduces local subdomain boundaries.
Definition communication.hpp:672
void pack_send_and_recv_local_subdomain_boundaries(const grid::shell::DistributedDomain &domain, const GridDataType &data, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_send_buffers, SubdomainNeighborhoodSendRecvBuffer< typename GridDataType::value_type, grid::grid_data_vec_dim< GridDataType >() > &boundary_recv_buffers)
Packs, sends and recvs local subdomain boundaries using two sets of buffers.
Definition communication.hpp:242
Definition boundary_mass.hpp:14
ScalarT supg_tau(const ScalarT vel_norm, const ScalarT kappa, const ScalarT h, const ScalarT Pe_tol=1e-8)
Definition unsteady_advection_diffusion_supg.hpp:18
constexpr void quad_felippa_3x2_quad_weights(T(&quad_weights)[quad_felippa_3x2_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:93
constexpr int quad_felippa_3x2_num_quad_points
Definition wedge/quadrature/quadrature.hpp:66
constexpr void quad_felippa_3x2_quad_points(dense::Vec< T, 3 >(&quad_points)[quad_felippa_3x2_num_quad_points])
Definition wedge/quadrature/quadrature.hpp:70
constexpr int num_nodes_per_wedge_surface
Definition kernel_helpers.hpp:6
void atomically_add_local_wedge_scalar_coefficients(const grid::Grid4DDataScalar< T > &global_coefficients, const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const dense::Vec< T, 6 >(&local_coefficients)[2])
Performs an atomic add of the two local wedge coefficient vectors of a hex cell into the global coeff...
Definition kernel_helpers.hpp:407
void wedge_surface_physical_coords(dense::Vec< T, 3 >(&wedge_surf_phy_coords)[num_wedges_per_hex_cell][num_nodes_per_wedge_surface], const grid::Grid3DDataVec< T, 3 > &lateral_grid, const int local_subdomain_id, const int x_cell, const int y_cell)
Extracts the (unit sphere) surface vertex coords of the two wedges of a hex cell.
Definition kernel_helpers.hpp:26
constexpr int num_wedges_per_hex_cell
Definition kernel_helpers.hpp:5
void extract_local_wedge_scalar_coefficients(dense::Vec< T, 6 >(&local_coefficients)[2], const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const grid::Grid4DDataScalar< T > &global_coefficients)
Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient v...
Definition kernel_helpers.hpp:306
void extract_local_wedge_vector_coefficients(dense::Vec< T, 6 >(&local_coefficients)[2], const int local_subdomain_id, const int x_cell, const int y_cell, const int r_cell, const int d, const grid::Grid4DDataVec< T, VecDim > &global_coefficients)
Extracts the local vector coefficients for the two wedges of a hex cell from the global coefficient v...
Definition kernel_helpers.hpp:356
constexpr int num_nodes_per_wedge
Definition kernel_helpers.hpp:7
constexpr dense::Vec< T, 3 > grad_shape(const int node_idx, const T xi, const T eta, const T zeta)
Gradient of the full shape function:
Definition integrands.hpp:228
constexpr T shape(const int node_idx, const T xi, const T eta, const T zeta)
(Tensor-product) Shape function.
Definition integrands.hpp:146
constexpr dense::Mat< T, 3, 3 > jac(const dense::Vec< T, 3 > &p1_phy, const dense::Vec< T, 3 > &p2_phy, const dense::Vec< T, 3 > &p3_phy, const T r_1, const T r_2, const T xi, const T eta, const T zeta)
Definition integrands.hpp:643
Kokkos::MDRangePolicy< Kokkos::Rank< 4 > > local_domain_md_range_policy_cells(const DistributedDomain &distributed_domain)
Definition spherical_shell.hpp:2668
Kokkos::View< ScalarType ***[VecDim], Layout > Grid3DDataVec
Definition grid_types.hpp:40
Kokkos::View< ScalarType ****[VecDim], Layout > Grid4DDataVec
Definition grid_types.hpp:43
Kokkos::View< ScalarType ****, Layout > Grid4DDataScalar
Definition grid_types.hpp:25
Kokkos::View< ScalarType **, Layout > Grid2DDataScalar
Definition grid_types.hpp:19
OperatorApplyMode
Modes for applying an operator to a vector.
Definition operator.hpp:30
@ Replace
Overwrite the destination vector.
OperatorCommunicationMode
Modes for communication during operator application.
Definition operator.hpp:40
@ CommunicateAdditively
Communicate and add results.
constexpr bool has_flag(E mask_value, E flag) noexcept
Checks if a bitmask value contains a specific flag.
Definition bit_masking.hpp:43
void fill(const T value)
Definition mat.hpp:201
constexpr Mat diagonal() const
Definition mat.hpp:377
static constexpr Mat diagonal_from_vec(const Vec< T, Rows > &diagonal)
Definition mat.hpp:79
Mat & hadamard_product(const Mat &mat)
Definition mat.hpp:213
T norm() const
Definition vec.hpp:76
void fill(const T value)
Definition vec.hpp:30