Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
euler_res_cpu.f90
Go to the documentation of this file.
1! Copyright (c) 2025-2026, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
38 use euler_residual, only : euler_rhs_t
39 use field, only : field_t
40 use ax_product, only : ax_t
41 use coefs, only : coef_t
42 use gather_scatter, only : gs_t
43 use num_types, only : rp
44 use operators, only : div, rotate_cyc
45 use gs_ops, only : gs_op_add
48 use field_list, only : field_list_t
49 implicit none
50 private
51
52 type, public, extends(euler_rhs_t) :: euler_res_cpu_t
53 contains
54 procedure, nopass :: step => advance_primitive_variables_cpu
55 procedure, nopass :: evaluate_rhs_cpu
56 end type euler_res_cpu_t
57
58contains
73 subroutine advance_primitive_variables_cpu(rho_field, m_x, m_y, m_z, &
74 E, p, u, v, w, Ax, &
75 coef, gs, h, effective_visc, rk_scheme, dt)
76 type(field_t), intent(inout) :: rho_field, m_x, m_y, m_z, E
77 type(field_t), intent(in) :: p, u, v, w, h, effective_visc
78 class(ax_t), intent(inout) :: Ax
79 type(coef_t), intent(inout) :: coef
80 type(gs_t), intent(inout) :: gs
81 class(runge_kutta_time_scheme_t), intent(in) :: rk_scheme
82 real(kind=rp), intent(in) :: dt
83 integer :: n, s, i, j, l
84 type(field_t), pointer :: k_rho_1, k_rho_2, k_rho_3, k_rho_4, &
85 k_m_x_1, k_m_x_2, k_m_x_3, k_m_x_4, &
86 k_m_y_1, k_m_y_2, k_m_y_3, k_m_y_4, &
87 k_m_z_1, k_m_z_2, k_m_z_3, k_m_z_4, &
88 k_E_1, k_E_2, k_E_3, k_E_4, &
89 temp_rho, temp_m_x, temp_m_y, temp_m_z, temp_E
90 integer :: tmp_indices(25)
91 type(field_list_t) :: k_rho, k_m_x, k_m_y, k_m_z, k_E
92 ! These contiguous pointers are necessary to ensure that
93 ! the Fujitsu compiler properly vectorizes the loop nests
94 real(kind=rp), contiguous, pointer :: k_rho_ptr(:,:,:,:)
95 real(kind=rp), contiguous, pointer :: k_m_x_ptr(:,:,:,:)
96 real(kind=rp), contiguous, pointer :: k_m_y_ptr(:,:,:,:)
97 real(kind=rp), contiguous, pointer :: k_m_z_ptr(:,:,:,:)
98 real(kind=rp), contiguous, pointer :: k_e_ptr(:,:,:,:)
99
100 n = p%dof%size()
101 s = rk_scheme%order
102 call neko_scratch_registry%request_field(k_rho_1, tmp_indices(1), .true.)
103 call neko_scratch_registry%request_field(k_rho_2, tmp_indices(2), .true.)
104 call neko_scratch_registry%request_field(k_rho_3, tmp_indices(3), .true.)
105 call neko_scratch_registry%request_field(k_rho_4, tmp_indices(4), .true.)
106 call neko_scratch_registry%request_field(k_m_x_1, tmp_indices(5), .true.)
107 call neko_scratch_registry%request_field(k_m_x_2, tmp_indices(6), .true.)
108 call neko_scratch_registry%request_field(k_m_x_3, tmp_indices(7), .true.)
109 call neko_scratch_registry%request_field(k_m_x_4, tmp_indices(8), .true.)
110 call neko_scratch_registry%request_field(k_m_y_1, tmp_indices(9), .true.)
111 call neko_scratch_registry%request_field(k_m_y_2, tmp_indices(10), .true.)
112 call neko_scratch_registry%request_field(k_m_y_3, tmp_indices(11), .true.)
113 call neko_scratch_registry%request_field(k_m_y_4, tmp_indices(12), .true.)
114 call neko_scratch_registry%request_field(k_m_z_1, tmp_indices(13), .true.)
115 call neko_scratch_registry%request_field(k_m_z_2, tmp_indices(14), .true.)
116 call neko_scratch_registry%request_field(k_m_z_3, tmp_indices(15), .true.)
117 call neko_scratch_registry%request_field(k_m_z_4, tmp_indices(16), .true.)
118 call neko_scratch_registry%request_field(k_e_1, tmp_indices(17), .true.)
119 call neko_scratch_registry%request_field(k_e_2, tmp_indices(18), .true.)
120 call neko_scratch_registry%request_field(k_e_3, tmp_indices(19), .true.)
121 call neko_scratch_registry%request_field(k_e_4, tmp_indices(20), .true.)
122 call neko_scratch_registry%request_field(temp_rho, tmp_indices(21), .false.)
123 call neko_scratch_registry%request_field(temp_m_x, tmp_indices(22), .false.)
124 call neko_scratch_registry%request_field(temp_m_y, tmp_indices(23), .false.)
125 call neko_scratch_registry%request_field(temp_m_z, tmp_indices(24), .false.)
126 call neko_scratch_registry%request_field(temp_e, tmp_indices(25), .false.)
127
128 ! Initialize Runge-Kutta stage variables for each conserved quantity
129 call k_rho%init(4)
130 call k_rho%assign(1, k_rho_1)
131 call k_rho%assign(2, k_rho_2)
132 call k_rho%assign(3, k_rho_3)
133 call k_rho%assign(4, k_rho_4)
134 call k_m_x%init(4)
135 call k_m_x%assign(1, k_m_x_1)
136 call k_m_x%assign(2, k_m_x_2)
137 call k_m_x%assign(3, k_m_x_3)
138 call k_m_x%assign(4, k_m_x_4)
139 call k_m_y%init(4)
140 call k_m_y%assign(1, k_m_y_1)
141 call k_m_y%assign(2, k_m_y_2)
142 call k_m_y%assign(3, k_m_y_3)
143 call k_m_y%assign(4, k_m_y_4)
144 call k_m_z%init(4)
145 call k_m_z%assign(1, k_m_z_1)
146 call k_m_z%assign(2, k_m_z_2)
147 call k_m_z%assign(3, k_m_z_3)
148 call k_m_z%assign(4, k_m_z_4)
149 call k_e%init(4)
150 call k_e%assign(1, k_e_1)
151 call k_e%assign(2, k_e_2)
152 call k_e%assign(3, k_e_3)
153 call k_e%assign(4, k_e_4)
154
155 ! Loop over Runge-Kutta stages. One parallel region per stage covers
156 ! both the initial copy and all (i-1) accumulation sweeps.
157 do i = 1, s
158 !$omp parallel private(j, l, k_rho_ptr) &
159 !$omp& private(k_m_x_ptr, k_m_y_ptr, k_m_z_ptr, k_E_ptr)
160
161 ! Copy current solution state to temporary arrays for this RK stage
162 !OCL NORECURRENCE, NOVREC, NOALIAS
163 !DIR$ CONCURRENT
164 !GCC$ ivdep
165 !$omp do simd
166 do l = 1, n
167 temp_rho%x(l,1,1,1) = rho_field%x(l,1,1,1)
168 temp_m_x%x(l,1,1,1) = m_x%x(l,1,1,1)
169 temp_m_y%x(l,1,1,1) = m_y%x(l,1,1,1)
170 temp_m_z%x(l,1,1,1) = m_z%x(l,1,1,1)
171 temp_e%x(l,1,1,1) = e%x(l,1,1,1)
172 end do
173 !$omp end do simd
174
175 ! Accumulate previous stage contributions using RK coefficients.
176 ! Each thread independently sets its private k_*_ptr aliases before
177 ! the worksharing loop; the implicit barrier at end of !$omp do simd
178 ! synchronises threads between j iterations.
179 do j = 1, i-1
180 k_rho_ptr => k_rho%items(j)%ptr%x
181 k_m_x_ptr => k_m_x%items(j)%ptr%x
182 k_m_y_ptr => k_m_y%items(j)%ptr%x
183 k_m_z_ptr => k_m_z%items(j)%ptr%x
184 k_e_ptr => k_e%items(j)%ptr%x
185
186 !OCL NORECURRENCE, NOVREC, NOALIAS
187 !DIR$ CONCURRENT
188 !GCC$ ivdep
189 !$omp do simd
190 do l = 1, n
191 temp_rho%x(l,1,1,1) = temp_rho%x(l,1,1,1) + &
192 dt * rk_scheme%coeffs_A(i, j) * k_rho_ptr(l,1,1,1)
193 temp_m_x%x(l,1,1,1) = temp_m_x%x(l,1,1,1) + &
194 dt * rk_scheme%coeffs_A(i, j) * k_m_x_ptr(l,1,1,1)
195 temp_m_y%x(l,1,1,1) = temp_m_y%x(l,1,1,1) + &
196 dt * rk_scheme%coeffs_A(i, j) * k_m_y_ptr(l,1,1,1)
197 temp_m_z%x(l,1,1,1) = temp_m_z%x(l,1,1,1) + &
198 dt * rk_scheme%coeffs_A(i, j) * k_m_z_ptr(l,1,1,1)
199 temp_e%x(l,1,1,1) = temp_e%x(l,1,1,1) + &
200 dt * rk_scheme%coeffs_A(i, j) * k_e_ptr(l,1,1,1)
201 end do
202 !$omp end do simd
203 end do
204 !$omp end parallel
205
206 ! Evaluate RHS terms for current stage using intermediate solution values
207 call evaluate_rhs_cpu(k_rho%items(i)%ptr, k_m_x%items(i)%ptr, &
208 k_m_y%items(i)%ptr, k_m_z%items(i)%ptr, &
209 k_e%items(i)%ptr, &
210 temp_rho, temp_m_x, temp_m_y, temp_m_z, temp_e, &
211 p, u, v, w, ax, &
212 coef, gs, h, effective_visc)
213 end do
214
215 ! Update the solution. Single parallel region covers all s stages.
216 !$omp parallel default(shared) private(i, l, k_rho_ptr) &
217 !$omp& private(k_m_x_ptr, k_m_y_ptr, k_m_z_ptr, k_E_ptr)
218 do i = 1, s
219 k_rho_ptr => k_rho%items(i)%ptr%x
220 k_m_x_ptr => k_m_x%items(i)%ptr%x
221 k_m_y_ptr => k_m_y%items(i)%ptr%x
222 k_m_z_ptr => k_m_z%items(i)%ptr%x
223 k_e_ptr => k_e%items(i)%ptr%x
224
225 !OCL NORECURRENCE, NOVREC, NOALIAS
226 !DIR$ CONCURRENT
227 !GCC$ ivdep
228 !$omp do simd
229 do l = 1, n
230 rho_field%x(l,1,1,1) = rho_field%x(l,1,1,1) + &
231 dt * rk_scheme%coeffs_b(i) * k_rho_ptr(l,1,1,1)
232 m_x%x(l,1,1,1) = m_x%x(l,1,1,1) + &
233 dt * rk_scheme%coeffs_b(i) * k_m_x_ptr(l,1,1,1)
234 m_y%x(l,1,1,1) = m_y%x(l,1,1,1) + &
235 dt * rk_scheme%coeffs_b(i) * k_m_y_ptr(l,1,1,1)
236 m_z%x(l,1,1,1) = m_z%x(l,1,1,1) + &
237 dt * rk_scheme%coeffs_b(i) * k_m_z_ptr(l,1,1,1)
238 e%x(l,1,1,1) = e%x(l,1,1,1) + &
239 dt * rk_scheme%coeffs_b(i) * k_e_ptr(l,1,1,1)
240 end do
241 !$omp end do simd
242 end do
243 !$omp end parallel
244
245 call neko_scratch_registry%relinquish_field(tmp_indices)
246
248
270 subroutine evaluate_rhs_cpu(rhs_rho_field, rhs_m_x, rhs_m_y, rhs_m_z, rhs_E, &
271 rho_field, m_x, m_y, m_z, E, p, u, v, w, Ax, &
272 coef, gs, h, effective_visc)
273 type(field_t), intent(inout) :: rhs_rho_field, &
274 rhs_m_x, rhs_m_y, rhs_m_z, rhs_e
275 type(field_t), intent(inout) :: rho_field, m_x, m_y, m_z, E
276 type(field_t), intent(in) :: p, u, v, w, h, effective_visc
277 class(ax_t), intent(inout) :: Ax
278 type(coef_t), intent(inout) :: coef
279 type(gs_t), intent(inout) :: gs
280 integer :: i, n
281 type(field_t), pointer :: f_x, f_y, f_z, &
282 visc_rho, visc_m_x, visc_m_y, visc_m_z, visc_E
283 integer :: tmp_indices(8)
284
285 n = coef%dof%size()
286 call neko_scratch_registry%request_field(f_x, tmp_indices(1), .false.)
287 call neko_scratch_registry%request_field(f_y, tmp_indices(2), .false.)
288 call neko_scratch_registry%request_field(f_z, tmp_indices(3), .false.)
289
290 ! Hoisted from below so the mult and h1=effective_visc sweeps can be
291 ! fused into one !$omp parallel do simd (registry calls are not
292 ! thread-safe so they cannot sit between two !$omp do constructs).
293 call neko_scratch_registry%request_field(visc_rho, tmp_indices(4), .false.)
294 call neko_scratch_registry%request_field(visc_m_x, tmp_indices(5), .false.)
295 call neko_scratch_registry%request_field(visc_m_y, tmp_indices(6), .false.)
296 call neko_scratch_registry%request_field(visc_m_z, tmp_indices(7), .false.)
297 call neko_scratch_registry%request_field(visc_e, tmp_indices(8), .false.)
298
300 ! Compute density flux divergence
301 call div(rhs_rho_field%x, m_x%x, m_y%x, m_z%x, coef)
302
304 ! Compute momentum flux divergences
305 ! m_x
306 !OCL NORECURRENCE, NOVREC, NOALIAS
307 !DIR$ CONCURRENT
308 !GCC$ ivdep
309 !$omp parallel do simd
310 do i = 1, n
311 f_x%x(i,1,1,1) = m_x%x(i,1,1,1) * m_x%x(i,1,1,1) / &
312 rho_field%x(i,1,1,1) + p%x(i,1,1,1)
313 f_y%x(i,1,1,1) = m_x%x(i,1,1,1) * m_y%x(i,1,1,1) / &
314 rho_field%x(i,1,1,1)
315 f_z%x(i,1,1,1) = m_x%x(i,1,1,1) * m_z%x(i,1,1,1) / &
316 rho_field%x(i,1,1,1)
317 end do
318 !$omp end parallel do simd
319 call div(rhs_m_x%x, f_x%x, f_y%x, f_z%x, coef)
320 ! m_y
321 !OCL NORECURRENCE, NOVREC, NOALIAS
322 !DIR$ CONCURRENT
323 !GCC$ ivdep
324 !$omp parallel do simd
325 do i = 1, n
326 f_x%x(i,1,1,1) = m_y%x(i,1,1,1) * m_x%x(i,1,1,1) / &
327 rho_field%x(i,1,1,1)
328 f_y%x(i,1,1,1) = m_y%x(i,1,1,1) * m_y%x(i,1,1,1) / &
329 rho_field%x(i,1,1,1) + p%x(i,1,1,1)
330 f_z%x(i,1,1,1) = m_y%x(i,1,1,1) * m_z%x(i,1,1,1) / &
331 rho_field%x(i,1,1,1)
332 end do
333 !$omp end parallel do simd
334 call div(rhs_m_y%x, f_x%x, f_y%x, f_z%x, coef)
335 ! m_z
336 !OCL NORECURRENCE, NOVREC, NOALIAS
337 !DIR$ CONCURRENT
338 !GCC$ ivdep
339 !$omp parallel do simd
340 do i = 1, n
341 f_x%x(i,1,1,1) = m_z%x(i,1,1,1) * m_x%x(i,1,1,1) / &
342 rho_field%x(i,1,1,1)
343 f_y%x(i,1,1,1) = m_z%x(i,1,1,1) * m_y%x(i,1,1,1) / &
344 rho_field%x(i,1,1,1)
345 f_z%x(i,1,1,1) = m_z%x(i,1,1,1) * m_z%x(i,1,1,1) / &
346 rho_field%x(i,1,1,1) + p%x(i,1,1,1)
347 end do
348 !$omp end parallel do simd
349 call div(rhs_m_z%x, f_x%x, f_y%x, f_z%x, coef)
350
352 ! Compute energy flux divergence
353 !OCL NORECURRENCE, NOVREC, NOALIAS
354 !DIR$ CONCURRENT
355 !GCC$ ivdep
356 !$omp parallel do simd
357 do i = 1, n
358 f_x%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * &
359 u%x(i,1,1,1)
360 f_y%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * &
361 v%x(i,1,1,1)
362 f_z%x(i,1,1,1) = (e%x(i,1,1,1) + p%x(i,1,1,1)) * &
363 w%x(i,1,1,1)
364 end do
365 !$omp end parallel do simd
366 call div(rhs_e%x, f_x%x, f_y%x, f_z%x, coef)
367
368 ! gs
369 call gs%op(rhs_rho_field, gs_op_add)
370 call rotate_cyc(rhs_m_x%x, rhs_m_y%x, rhs_m_z%x, 1, coef)
371 call gs%op(rhs_m_x, gs_op_add)
372 call gs%op(rhs_m_y, gs_op_add)
373 call gs%op(rhs_m_z, gs_op_add)
374 call rotate_cyc(rhs_m_x%x, rhs_m_y%x, rhs_m_z%x, 0, coef)
375 call gs%op(rhs_e, gs_op_add)
376
377 ! Apply multiplicity to the inviscid RHS and set h1 to the effective
378 ! viscosity for the Laplacian. Fused: one fork-join instead of two, and
379 ! coef%mult / effective_visc share L1.
380 !OCL NORECURRENCE, NOVREC, NOALIAS
381 !DIR$ CONCURRENT
382 !GCC$ ivdep
383 !$omp parallel do simd
384 do i = 1, n
385 rhs_rho_field%x(i,1,1,1) = rhs_rho_field%x(i,1,1,1) * &
386 coef%mult(i,1,1,1)
387 rhs_m_x%x(i,1,1,1) = rhs_m_x%x(i,1,1,1) * coef%mult(i,1,1,1)
388 rhs_m_y%x(i,1,1,1) = rhs_m_y%x(i,1,1,1) * coef%mult(i,1,1,1)
389 rhs_m_z%x(i,1,1,1) = rhs_m_z%x(i,1,1,1) * coef%mult(i,1,1,1)
390 rhs_e%x(i,1,1,1) = rhs_e%x(i,1,1,1) * coef%mult(i,1,1,1)
391 coef%h1(i,1,1,1) = effective_visc%x(i,1,1,1)
392 end do
393 !$omp end parallel do simd
394
395 ! Calculate artificial diffusion with variable viscosity
396 call ax%compute(visc_rho%x, rho_field%x, coef, p%msh, p%Xh)
397 call ax%compute(visc_m_x%x, m_x%x, coef, p%msh, p%Xh)
398 call ax%compute(visc_m_y%x, m_y%x, coef, p%msh, p%Xh)
399 call ax%compute(visc_m_z%x, m_z%x, coef, p%msh, p%Xh)
400 call ax%compute(visc_e%x, e%x, coef, p%msh, p%Xh)
401
402 ! gs. h1=1.0 reset is deferred into the final accumulation below; safe
403 ! because nothing here (gs%op, rotate_cyc) reads coef%h1.
404 call gs%op(visc_rho, gs_op_add)
405 call rotate_cyc(visc_m_x%x, visc_m_y%x, visc_m_z%x, 1, coef)
406 call gs%op(visc_m_x, gs_op_add)
407 call gs%op(visc_m_y, gs_op_add)
408 call gs%op(visc_m_z, gs_op_add)
409 call rotate_cyc(visc_m_x%x, visc_m_y%x, visc_m_z%x, 0, coef)
410 call gs%op(visc_e, gs_op_add)
411
412 ! Move div to the rhs, apply artificial viscosity, and reset h1 to 1.
413 ! The viscosity coefficient is already included in the Laplacian operator.
414 ! Fused: one fork-join instead of two, and coef%Binv stays in L1.
415 !OCL NORECURRENCE, NOVREC, NOALIAS
416 !DIR$ CONCURRENT
417 !GCC$ ivdep
418 !$omp parallel do simd
419 do i = 1, n
420 rhs_rho_field%x(i,1,1,1) = -rhs_rho_field%x(i,1,1,1) - &
421 coef%Binv(i,1,1,1) * visc_rho%x(i,1,1,1)
422 rhs_m_x%x(i,1,1,1) = -rhs_m_x%x(i,1,1,1) - &
423 coef%Binv(i,1,1,1) * visc_m_x%x(i,1,1,1)
424 rhs_m_y%x(i,1,1,1) = -rhs_m_y%x(i,1,1,1) - &
425 coef%Binv(i,1,1,1) * visc_m_y%x(i,1,1,1)
426 rhs_m_z%x(i,1,1,1) = -rhs_m_z%x(i,1,1,1) - &
427 coef%Binv(i,1,1,1) * visc_m_z%x(i,1,1,1)
428 rhs_e%x(i,1,1,1) = -rhs_e%x(i,1,1,1) - &
429 coef%Binv(i,1,1,1) * visc_e%x(i,1,1,1)
430 coef%h1(i,1,1,1) = 1.0_rp
431 end do
432 !$omp end parallel do simd
433
434 call neko_scratch_registry%relinquish_field(tmp_indices)
435 end subroutine evaluate_rhs_cpu
436
437end module euler_res_cpu
Compute the divergence of a vector field.
Definition operators.f90:85
Apply cyclic boundary condition to a vector field.
Defines a Matrix-vector product.
Definition ax.f90:34
Coefficients.
Definition coef.f90:34
This implements CPU-based residual calculations for the Euler equations. It handles the time advancem...
subroutine evaluate_rhs_cpu(rhs_rho_field, rhs_m_x, rhs_m_y, rhs_m_z, rhs_e, rho_field, m_x, m_y, m_z, e, p, u, v, w, ax, coef, gs, h, effective_visc)
Evaluates the right-hand side of the Euler equations Computes fluxes, divergence terms,...
subroutine advance_primitive_variables_cpu(rho_field, m_x, m_y, m_z, e, p, u, v, w, ax, coef, gs, h, effective_visc, rk_scheme, dt)
Advances the primitive variables (density, momentum, energy) in time using a Runge-Kutta scheme.
Defines a field.
Definition field.f90:34
Gather-scatter.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Base type for a matrix-vector product providing .
Definition ax.f90:43
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Abstract type to compute rhs.
Definition euler_res.f90:44
field_list_t, To be able to group fields together