Loading...
Searching...
No Matches
unsteady_advection_diffusion_supg.hpp
Go to the documentation of this file.
1
2#pragma once
3
4#include "../../quadrature/quadrature.hpp"
6#include "dense/vec.hpp"
10#include "linalg/operator.hpp"
11#include "linalg/vector.hpp"
12#include "linalg/vector_q1.hpp"
13
15
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 )
19{
20 const ScalarT kappa_min = ScalarT( 1e-8 );
21 const ScalarT kappa_virtual = Kokkos::max( kappa, kappa_min );
22
23 // Guarding against zero velocity.
24 const ScalarT eps_vel = ScalarT( 1e-12 );
25
26 if ( vel_norm <= eps_vel )
27 {
28 return 0.0;
29 }
30
31 // For small Peclet numbers (diffusion-dominated flow) set tau to 0 (no stabilization).
32 const ScalarT Pe = vel_norm * h / ( 2.0 * kappa_virtual );
33
34 if ( Pe <= Pe_tol )
35 {
36 return 0.0;
37 }
38
39 // Clamp tau to avoid huge values
40 const ScalarT tau_max = ScalarT( 10.0 ) * h / Kokkos::max( vel_norm, eps_vel );
41
42 // Finally, compute tau.
43 // Note: coth(Pe) = 1/tanh(Pe)
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;
46
47 return ( tau > tau_max ) ? tau_max : tau;
48 // return tau;
49}
50
51/// \brief Linear operator for a method-of-lines discretization of the unsteady advection-diffusion equation with SUPG
52/// (streamline upwind Petrov-Galerkin) stabilization.
53///
54/// # Continuous problem
55///
56/// The unsteady advection-diffusion equation is given by
57/// \f[
58/// \frac{\partial}{\partial t}T + \mathbf{u} \cdot \nabla T - \nabla \cdot (\kappa \nabla T) = f
59/// \f]
60/// where \f$ T \f$ is the scalar temperature solution, \f$ \mathbf{u} \f$ a given velocity field,
61/// \f$ f \f$ a given source term, and \f$ \kappa \f$ a given diffusivity function.
62///
63/// # Space discretization
64///
65/// \note We assume here that we have \f$ \kappa|_E = \mathrm{const} \f$ on each element \f$E\f$ which simplifies the
66/// implementation of the SUPG stabilization for linear finite elements as certain terms drop. Currently, \f$\kappa =
67/// \mathrm{const}\f$ globally, but once that is changed, we will need to average on every element it in the kernel or
68/// pass it in as an elementwise constant (FE) function.
69///
70/// We first discretize in space, then in time (method of lines). After discretization in space, we get the system of
71/// ODEs in time
72/// \f[
73/// M \frac{d}{dt}T + (C + K + G)T = F + F_\mathrm{SUPG}
74/// \f]
75/// where
76///
77/// Term | Description | Bilinear form
78/// -------------|--------------|--------------
79/// \f$M_{ij}\f$ | mass | \f$ \int \phi_i \phi_j \f$
80/// \f$C_{ij}\f$ | advection | \f$ \int \phi_i (\mathbf{u} \cdot \nabla \phi_j) \f$
81/// \f$K_{ij}\f$ | diffusion | \f$ \int \nabla \phi_i \cdot (\kappa \nabla \phi_j) \f$
82/// \f$G_{ij}\f$ | SUPG adv-adv | \f$ \sum_E \int_E \tau_E (\mathbf{u} \cdot \nabla \phi_i) (\mathbf{u} \cdot \nabla \phi_j) \f$
83/// \f$F_i\f$ | forcing | \f$ \int \phi_i f \f$
84/// \f$(F_\mathrm{SUPG})_{i}\f$ | SUPG forcing | \f$ \sum_E \int_E \tau_E (\mathbf{u} \cdot \nabla \phi_i) f \f$.
85///
86/// \note After brief testing, it seems that in general the term \f$F_\mathrm{SUPG}\f$ does not always improve the
87/// computed solutions. It even appears to slightly increase the error sometimes. So setting
88/// \f$F_\mathrm{SUPG} = 0\f$ can really be just fine in certain settings.
89///
90/// # Time discretization
91///
92/// For the time discretization, we employ implicit BDF schemes. A general formula is
93///
94/// \f[
95/// (\alpha_0 M + \Delta t A) T^{n+1} = - \sum_{j=1}^k \alpha_j M T^{n+1-j} + \Delta t R^{n+1}
96/// \f]
97///
98/// where \f$A = (C + K + G)\f$ and \f$R = F + F_{\mathrm{SUPG}}\f$.
99///
100/// We recover the common BDF1 (backward or implicit Euler) and BDF2 schemes by choosing:
101///
102/// Scheme | \f$k\f$ | \f$\alpha\f$ | full equation
103/// --------------------------------|---------|--------------------------------------|--------------
104/// BDF1 (backward/implicit Euler) | 1 | \f$[1, -1]\f$ | \f$(M + \Delta t A) T^{n+1} = M T^{n} + \Delta t R^{n+1}\f$
105/// BDF2 | 2 | \f$[\frac{3}{2}, -2, \frac{1}{2}]\f$ | \f$(\frac{3}{2} M + \Delta t A) T^{n+1} = 2 M T^{n} - \frac{1}{2} M T^{n-1} + \Delta t R^{n+1}\f$
106///
107/// The RHS term must be built with appropriate other classes/function. This class is only concerned with the
108/// matrix-free evaluation of the LHS system matrix
109/// \f[
110/// \alpha_0 M + \Delta t A.
111/// \f]
112///
113/// The parameters \f$\alpha_0\f$ and \f$\Delta t\f$ are set through the constructor via `mass_scaling` and `dt`
114/// respectively.
115///
116/// # SUPG stabilization
117///
118/// Several choices for the stabilization parameter \f$ \tau_E \f$ are possible.
119/// As it is commonly selected in the literature (see, e.g., John and Knobloch (2006)), we set on element \f$E\f$
120/// \f[
121/// \tau_E = \frac{h_E}{2 \|\mathbf{u}\|} \left(\coth(\mathrm{Pe_E}) - \frac{1}{\mathrm{Pe_E}}\right)
122/// \f]
123/// where
124/// \f[
125/// \mathrm{Pe_E} = \frac{\|\mathbf{u}\|h_E}{2 \kappa}
126/// \f]
127/// with some precautionary measures to avoid edge cases that would result in division by zero etc.
128///
129template < typename ScalarT, int VelocityVecDim = 3 >
131{
132 public:
135 using ScalarType = ScalarT;
136
137 private:
139
143
145
146 ScalarT diffusivity_;
147 ScalarT dt_;
148
149 bool treat_boundary_;
150 bool diagonal_;
151 ScalarT mass_scaling_;
152 bool lumped_mass_;
153
154 linalg::OperatorApplyMode operator_apply_mode_;
155 linalg::OperatorCommunicationMode operator_communication_mode_;
156
159
163
164 public:
166 const grid::shell::DistributedDomain& domain,
171 const ScalarT diffusivity,
172 const ScalarT dt,
173 bool treat_boundary,
174 bool diagonal = false,
175 ScalarT mass_scaling = 1.0,
176 bool lumped_mass = false,
178 linalg::OperatorCommunicationMode operator_communication_mode =
180 : domain_( domain )
181 , grid_( grid )
182 , radii_( radii )
183 , boundary_mask_( boundary_mask )
184 , velocity_( velocity )
185 , diffusivity_( diffusivity )
186 , dt_( dt )
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 )
193 // TODO: we can reuse the send and recv buffers and pass in from the outside somehow
194 , send_buffers_( domain )
195 , recv_buffers_( domain )
196 {}
197
198 ScalarT& dt() { return dt_; }
199 const ScalarT& dt() const { return dt_; }
200
201 void apply_impl( const SrcVectorType& src, DstVectorType& dst )
202 {
203 util::Timer timer_apply( "ad_supg_apply" );
204
205 if ( operator_apply_mode_ == linalg::OperatorApplyMode::Replace )
206 {
207 assign( dst, 0 );
208 }
209
210 src_ = src.grid_data();
211 dst_ = dst.grid_data();
212 vel_grid_ = velocity_.grid_data();
213
214 util::Timer timer_kernel( "ad_supg_kernel" );
215 Kokkos::parallel_for( "matvec", grid::shell::local_domain_md_range_policy_cells( domain_ ), *this );
216 Kokkos::fence();
217 timer_kernel.stop();
218
219 if ( operator_communication_mode_ == linalg::OperatorCommunicationMode::CommunicateAdditively )
220 {
221 util::Timer timer_comm( "ad_supg_comm" );
222
224 domain_, dst_, send_buffers_, recv_buffers_ );
226 }
227 }
228
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
231 {
232 // Gather surface points for each wedge.
234 wedge_surface_physical_coords( wedge_phy_surf, grid_, local_subdomain_id, x_cell, y_cell );
235
236 // Gather wedge radii.
237 const ScalarT r_1 = radii_( local_subdomain_id, r_cell );
238 const ScalarT r_2 = radii_( local_subdomain_id, r_cell + 1 );
239
240 // Quadrature points.
241 constexpr auto num_quad_points = quadrature::quad_felippa_3x2_num_quad_points;
242
243 dense::Vec< ScalarT, 3 > quad_points[num_quad_points];
244 ScalarT quad_weights[num_quad_points];
245
248
249 // Interpolating velocity into quad points.
250
252 dense::Vec< ScalarT, 6 > vel_coeffs[VelocityVecDim][num_wedges_per_hex_cell];
253
254 for ( int d = 0; d < VelocityVecDim; d++ )
255 {
257 vel_coeffs[d], local_subdomain_id, x_cell, y_cell, r_cell, d, vel_grid_ );
258 }
259
260 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
261 {
262 for ( int q = 0; q < num_quad_points; q++ )
263 {
264 for ( int i = 0; i < num_nodes_per_wedge; i++ )
265 {
266 const auto shape_i = shape( i, quad_points[q] );
267 for ( int d = 0; d < VelocityVecDim; d++ )
268 {
269 vel_interp[wedge][q]( d ) += vel_coeffs[d][wedge]( i ) * shape_i;
270 }
271 }
272 }
273 }
274
275 // Let's compute the streamline diffusivity.
276
277 ScalarT streamline_diffusivity[num_wedges_per_hex_cell];
278
279 // Far from accurate but for now assume h = r.
280 const auto h = r_2 - r_1;
281
282 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
283 {
284 ScalarT tau_accum = 0.0;
285 ScalarT waccum = 0.0;
286
287 for ( int q = 0; q < num_quad_points; ++q )
288 {
289 // get velocity at this quad point
290 const auto& uq = vel_interp[wedge][q];
291 const ScalarT vel_norm_q = uq.norm();
292
293 const ScalarT tau_q = supg_tau< ScalarT >( vel_norm_q, diffusivity_, h, 1e-08 );
294
295 // quadrature weight for this point (if you have weights)
296 const ScalarT wq = quad_weights[q]; // if not available, use 1.0
297 tau_accum += tau_q * wq;
298 waccum += wq;
299 }
300
301 // final cell/wedge tau: volume-weighted average
302 ScalarT tau_cell = ( waccum > 0.0 ) ? ( tau_accum / waccum ) : 0.0;
303 streamline_diffusivity[wedge] = tau_cell;
304 }
305
306 // Compute the local element matrix.
309
310 for ( int q = 0; q < num_quad_points; q++ )
311 {
312 const auto w = quad_weights[q];
313
314 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
315 {
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();
319
320 const auto vel = vel_interp[wedge][q];
321
322 for ( int i = 0; i < num_nodes_per_wedge; i++ )
323 {
324 const auto shape_i = shape( i, quad_points[q] );
325 const auto grad_i = J_inv_transposed * grad_shape( i, quad_points[q] );
326
327 for ( int j = 0; j < num_nodes_per_wedge; j++ )
328 {
329 const auto shape_j = shape( j, quad_points[q] );
330 const auto grad_j = J_inv_transposed * grad_shape( j, quad_points[q] );
331
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 ) );
337
338 M[wedge]( i, j ) += w * mass_scaling_ * mass * det;
339 A[wedge]( i, j ) += w * dt_ * ( diffusion + advection + streamline_diffusion ) * det;
340 }
341 }
342 }
343 }
344
345 if ( lumped_mass_ )
346 {
348 ones.fill( 1.0 );
351 }
352
353 if ( treat_boundary_ )
354 {
355 const int at_cmb_boundary = util::has_flag(
356 boundary_mask_( local_subdomain_id, x_cell, y_cell, r_cell ), grid::shell::ShellBoundaryFlag::CMB );
357 const int at_surface_boundary = util::has_flag(
358 boundary_mask_( local_subdomain_id, x_cell, y_cell, r_cell + 1 ),
360
361 for ( int wedge = 0; wedge < num_wedges_per_hex_cell; wedge++ )
362 {
363 dense::Mat< ScalarT, 6, 6 > boundary_mask;
364 boundary_mask.fill( 1.0 );
365 if ( at_cmb_boundary )
366 {
367 // Inner boundary (CMB).
368 for ( int i = 0; i < 6; i++ )
369 {
370 for ( int j = 0; j < 6; j++ )
371 {
372 if ( i != j && ( i < 3 || j < 3 ) )
373 {
374 boundary_mask( i, j ) = 0.0;
375 }
376 }
377 }
378 }
379
380 if ( at_surface_boundary )
381 {
382 // Outer boundary (surface).
383 for ( int i = 0; i < 6; i++ )
384 {
385 for ( int j = 0; j < 6; j++ )
386 {
387 if ( i != j && ( i >= 3 || j >= 3 ) )
388 {
389 boundary_mask( i, j ) = 0.0;
390 }
391 }
392 }
393 }
394
395 M[wedge].hadamard_product( boundary_mask );
396 A[wedge].hadamard_product( boundary_mask );
397 }
398 }
399
400 if ( diagonal_ )
401 {
402 M[0] = M[0].diagonal();
403 M[1] = M[1].diagonal();
404 A[0] = A[0].diagonal();
405 A[1] = A[1].diagonal();
406 }
407
409 extract_local_wedge_scalar_coefficients( src, local_subdomain_id, x_cell, y_cell, r_cell, src_ );
410
412
413 dst[0] = ( M[0] + A[0] ) * src[0];
414 dst[1] = ( M[1] + A[1] ) * src[1];
415
416 atomically_add_local_wedge_scalar_coefficients( dst_, local_subdomain_id, x_cell, y_cell, r_cell, dst );
417 }
418};
419
421
422} // namespace terra::fe::wedge::operators::shell
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
Definition mat.hpp:10
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