Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_volflow.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2! Copyright (c) 2025-2026, The Neko Authors
3!
4! The UChicago Argonne, LLC as Operator of Argonne National
5! Laboratory holds copyright in the Software. The copyright holder
6! reserves all rights except those expressly granted to licensees,
7! and U.S. Government license rights.
8!
9! Redistribution and use in source and binary forms, with or without
10! modification, are permitted provided that the following conditions
11! are met:
12!
13! 1. Redistributions of source code must retain the above copyright
14! notice, this list of conditions and the disclaimer below.
15!
16! 2. Redistributions in binary form must reproduce the above copyright
17! notice, this list of conditions and the disclaimer (as noted below)
18! in the documentation and/or other materials provided with the
19! distribution.
20!
21! 3. Neither the name of ANL nor the names of its contributors
22! may be used to endorse or promote products derived from this software
23! without specific prior written permission.
24!
25! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
26! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
27! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
28! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
29! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
30! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
31! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
32! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
33! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
34! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
35! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
36! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37!
38! Additional BSD Notice
39! ---------------------
40! 1. This notice is required to be provided under our contract with
41! the U.S. Department of Energy (DOE). This work was produced at
42! Argonne National Laboratory under Contract
43! No. DE-AC02-06CH11357 with the DOE.
44!
45! 2. Neither the United States Government nor UCHICAGO ARGONNE,
46! LLC nor any of their employees, makes any warranty,
47! express or implied, or assumes any liability or responsibility for the
48! accuracy, completeness, or usefulness of any information, apparatus,
49! product, or process disclosed, or represents that its use would not
50! infringe privately-owned rights.
51!
52! 3. Also, reference herein to any specific commercial products, process,
53! or services by trade name, trademark, manufacturer or otherwise does
54! not necessarily constitute or imply its endorsement, recommendation,
55! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
56! The views and opinions of authors expressed
57! herein do not necessarily state or reflect those of the United States
58! Government or UCHICAGO ARGONNE, LLC, and shall
59! not be used for advertising or product endorsement purposes.
60!
62 use operators, only : opgrad, cdtp, rotate_cyc
63 use num_types, only : rp
64 use mathops, only : opchsign
65 use krylov, only : ksp_t, ksp_monitor_t
66 use precon, only : pc_t
67 use dofmap, only : dofmap_t
68 use field, only : field_t
69 use coefs, only : coef_t
70 use time_state, only : time_state_t
72 use math, only : copy, glsc2, glmin, glmax, add2, abscmp
77 use gather_scatter, only : gs_t, gs_op_add
78 use json_module, only : json_file
81 use bc_list, only : bc_list_t
82 use ax_product, only : ax_t
84 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_sum
85 use logger, only : log_size, neko_log
86 implicit none
87 private
88
90 type, public :: fluid_volflow_t
91 integer :: flow_dir
92 logical :: avflow
93 logical :: log = .true.
94 real(kind=rp) :: flow_rate
95 real(kind=rp) :: dtlag = 0d0
96 real(kind=rp) :: bdlag = 0d0
97 type(field_t) :: u_vol, v_vol, w_vol, p_vol
98 real(kind=rp) :: domain_length, base_flow
99 contains
100 procedure, pass(this) :: init => fluid_vol_flow_init
101 procedure, pass(this) :: free => fluid_vol_flow_free
102 procedure, pass(this) :: adjust => fluid_vol_flow
103 procedure, private, pass(this) :: compute => fluid_vol_flow_compute
104 end type fluid_volflow_t
105
106contains
107
108 subroutine fluid_vol_flow_init(this, dm_Xh, params)
109 class(fluid_volflow_t), intent(inout) :: this
110 type(dofmap_t), target, intent(in) :: dm_Xh
111 type(json_file), intent(inout) :: params
112 logical average, log_output
113 integer :: direction
114 real(kind=rp) :: rate
115
116 call this%free()
117
118 !Initialize vol_flow (if there is a forced volume flow)
119 call json_get_or_lookup(params, 'case.fluid.flow_rate_force.direction', &
120 direction)
121 call json_get_or_lookup(params, 'case.fluid.flow_rate_force.value', rate)
122 call json_get(params, 'case.fluid.flow_rate_force.use_averaged_flow',&
123 average)
124 call json_get_or_default(params, 'case.fluid.flow_rate_force.log', &
125 log_output, .true.)
126
127 this%flow_dir = direction
128 this%avflow = average
129 this%log = log_output
130 this%flow_rate = rate
131
132 if (this%flow_dir .ne. 0) then
133 call this%u_vol%init(dm_xh, 'u_vol')
134 call this%v_vol%init(dm_xh, 'v_vol')
135 call this%w_vol%init(dm_xh, 'w_vol')
136 call this%p_vol%init(dm_xh, 'p_vol')
137 end if
138
139 end subroutine fluid_vol_flow_init
140
141 subroutine fluid_vol_flow_free(this)
142 class(fluid_volflow_t), intent(inout) :: this
143
144 call this%u_vol%free()
145 call this%v_vol%free()
146 call this%w_vol%free()
147 call this%p_vol%free()
148
149 end subroutine fluid_vol_flow_free
150
154 subroutine fluid_vol_flow_compute(this, u_res, v_res, w_res, p_res, &
155 ext_bdf, gs_Xh, c_Xh, rho, mu, bd, dt, &
156 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
157 Ax_vel, Ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
158 vel_max_iter)
159 class(fluid_volflow_t), intent(inout) :: this
160 type(field_t), intent(inout) :: u_res, v_res, w_res, p_res
161 type(coef_t), intent(inout) :: c_Xh
162 type(gs_t), intent(inout) :: gs_Xh
163 type(time_scheme_controller_t), intent(in) :: ext_bdf
164 type(bc_list_t), intent(inout) :: bclst_dp, bclst_du, bclst_dv, bclst_dw
165 type(bc_list_t), intent(inout) :: bclst_vel_res
166 class(ax_t), intent(in) :: Ax_vel
167 class(ax_t), intent(in) :: Ax_prs
168 class(ksp_t), intent(inout) :: ksp_prs, ksp_vel
169 class(pc_t), intent(inout) :: pc_prs, pc_vel
170 real(kind=rp), intent(in) :: bd
171 real(kind=rp), intent(in) :: rho, dt
172 type(field_t) :: mu
173 integer, intent(in) :: vel_max_iter, prs_max_iter
174 integer :: n, i
175 real(kind=rp) :: xlmin, xlmax
176 real(kind=rp) :: ylmin, ylmax
177 real(kind=rp) :: zlmin, zlmax
178 type(ksp_monitor_t) :: ksp_results(4)
179 type(field_t), pointer :: ta1, ta2, ta3
180 integer :: temp_indices(3)
181
182 call neko_scratch_registry%request_field(ta1, temp_indices(1), .false.)
183 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
184 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
185
186
187 associate(msh => c_xh%msh, p_vol => this%p_vol, &
188 u_vol => this%u_vol, v_vol => this%v_vol, w_vol => this%w_vol)
189
190 n = c_xh%dof%size()
191 xlmin = glmin(c_xh%dof%x, n)
192 xlmax = glmax(c_xh%dof%x, n)
193 ylmin = glmin(c_xh%dof%y, n) ! for Y!
194 ylmax = glmax(c_xh%dof%y, n)
195 zlmin = glmin(c_xh%dof%z, n) ! for Z!
196 zlmax = glmax(c_xh%dof%z, n)
197 if (this%flow_dir .eq. 1) then
198 this%domain_length = xlmax - xlmin
199 end if
200 if (this%flow_dir .eq. 2) then
201 this%domain_length = ylmax - ylmin
202 end if
203 if (this%flow_dir .eq. 3) then
204 this%domain_length = zlmax - zlmin
205 end if
206
207 if (neko_bcknd_device .eq. 1) then
208 call device_cfill(c_xh%h1_d, 1.0_rp/rho, n)
209 call device_rzero(c_xh%h2_d, n)
210 else
211 do i = 1, n
212 c_xh%h1(i,1,1,1) = 1.0_rp / rho
213 c_xh%h2(i,1,1,1) = 0.0_rp
214 end do
215 end if
216 c_xh%ifh2 = .false.
217
218 ! Compute pressure
219
220 if (this%flow_dir .eq. 1) then
221 call cdtp(p_res%x, c_xh%h1, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
222 end if
223
224 if (this%flow_dir .eq. 2) then
225 call cdtp(p_res%x, c_xh%h1, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
226 end if
227
228 if (this%flow_dir .eq. 3) then
229 call cdtp(p_res%x, c_xh%h1, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
230 end if
231
232 call gs_xh%op(p_res, gs_op_add)
233 call bclst_dp%apply_scalar(p_res%x, n)
234 call pc_prs%update()
235 ksp_results(1) = ksp_prs%solve(ax_prs, p_vol, p_res%x, n, &
236 c_xh, bclst_dp, gs_xh, prs_max_iter)
237
238 ! Compute velocity
239
240 call opgrad(u_res%x, v_res%x, w_res%x, p_vol%x, c_xh)
241
242 if (neko_bcknd_device .eq. 1) then
243 call device_opchsign(u_res%x_d, v_res%x_d, w_res%x_d, msh%gdim, n)
244 call device_copy(ta1%x_d, c_xh%B_d, n)
245 call device_copy(ta2%x_d, c_xh%B_d, n)
246 call device_copy(ta3%x_d, c_xh%B_d, n)
247 else
248 call opchsign(u_res%x, v_res%x, w_res%x, msh%gdim, n)
249 call copy(ta1%x, c_xh%B, n)
250 call copy(ta2%x, c_xh%B, n)
251 call copy(ta3%x, c_xh%B, n)
252 end if
253 call bclst_vel_res%apply_vector(ta1%x, ta2%x, ta3%x, n)
254
255 ! add forcing
256
257 if (neko_bcknd_device .eq. 1) then
258 if (this%flow_dir .eq. 1) then
259 call device_add2(u_res%x_d, ta1%x_d, n)
260 else if (this%flow_dir .eq. 2) then
261 call device_add2(v_res%x_d, ta2%x_d, n)
262 else if (this%flow_dir .eq. 3) then
263 call device_add2(w_res%x_d, ta3%x_d, n)
264 end if
265 else
266 if (this%flow_dir .eq. 1) then
267 call add2(u_res%x, ta1%x, n)
268 else if (this%flow_dir .eq. 2) then
269 call add2(v_res%x, ta2%x, n)
270 else if (this%flow_dir .eq. 3) then
271 call add2(w_res%x, ta3%x, n)
272 end if
273 end if
274
275 if (neko_bcknd_device .eq. 1) then
276 call device_copy(c_xh%h1_d, mu%x_d, n)
277 call device_cfill(c_xh%h2_d, rho * (bd / dt), n)
278 else
279 call copy(c_xh%h1, mu%x, n)
280 c_xh%h2 = rho * (bd / dt)
281 end if
282 c_xh%ifh2 = .true.
283
284 call rotate_cyc(u_res%x, v_res%x, w_res%x, 1, c_xh)
285 call gs_xh%op(u_res, gs_op_add)
286 call gs_xh%op(v_res, gs_op_add)
287 call gs_xh%op(w_res, gs_op_add)
288 call rotate_cyc(u_res%x, v_res%x, w_res%x, 0, c_xh)
289
290 call bclst_vel_res%apply_vector(u_res%x, v_res%x, w_res%x, n)
291 call pc_vel%update()
292
293 ksp_results(2:4) = ksp_vel%solve_coupled(ax_vel, &
294 u_vol, v_vol, w_vol, &
295 u_res%x, v_res%x, w_res%x, &
296 n, c_xh, &
297 bclst_du, bclst_dv, bclst_dw, &
298 gs_xh, vel_max_iter)
299
300 if (neko_bcknd_device .eq. 1) then
301 if (this%flow_dir .eq. 1) then
302 this%base_flow = &
303 device_glsc2(u_vol%x_d, c_xh%B_d, n) / this%domain_length
304 end if
305
306 if (this%flow_dir .eq. 2) then
307 this%base_flow = &
308 device_glsc2(v_vol%x_d, c_xh%B_d, n) / this%domain_length
309 end if
310
311 if (this%flow_dir .eq. 3) then
312 this%base_flow = &
313 device_glsc2(w_vol%x_d, c_xh%B_d, n) / this%domain_length
314 end if
315 else
316 if (this%flow_dir .eq. 1) then
317 this%base_flow = glsc2(u_vol%x, c_xh%B, n) / this%domain_length
318 end if
319
320 if (this%flow_dir .eq. 2) then
321 this%base_flow = glsc2(v_vol%x, c_xh%B, n) / this%domain_length
322 end if
323
324 if (this%flow_dir .eq. 3) then
325 this%base_flow = glsc2(w_vol%x, c_xh%B, n) / this%domain_length
326 end if
327 end if
328 end associate
329
330 call neko_scratch_registry%relinquish_field(temp_indices)
331 end subroutine fluid_vol_flow_compute
332
342 subroutine fluid_vol_flow(this, u, v, w, p, u_res, v_res, w_res, p_res, &
343 c_Xh, gs_Xh, ext_bdf, rho, mu, dt, time, &
344 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
345 Ax_vel, Ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
346 vel_max_iter)
347
348 class(fluid_volflow_t), intent(inout) :: this
349 type(field_t), intent(inout) :: u, v, w, p
350 type(field_t), intent(inout) :: u_res, v_res, w_res, p_res
351 type(coef_t), intent(inout) :: c_Xh
352 type(gs_t), intent(inout) :: gs_Xh
353 type(time_scheme_controller_t), intent(in) :: ext_bdf
354 type(time_state_t), intent(in) :: time
355 real(kind=rp), intent(in) :: rho, dt
356 type(field_t) :: mu
357 type(bc_list_t), intent(inout) :: bclst_dp, bclst_du, bclst_dv, bclst_dw
358 type(bc_list_t), intent(inout) :: bclst_vel_res
359 class(ax_t), intent(in) :: Ax_vel
360 class(ax_t), intent(in) :: Ax_prs
361 class(ksp_t), intent(inout) :: ksp_prs, ksp_vel
362 class(pc_t), intent(inout) :: pc_prs, pc_vel
363 integer, intent(in) :: prs_max_iter, vel_max_iter
364 real(kind=rp) :: ifcomp, flow_rate, xsec
365 real(kind=rp) :: current_flow, delta_flow, scale
366 integer :: n, ierr, i
367 character(len=5) :: flow_dir_label
368 character(len=12) :: step_str
369
370 character(len=200) :: log_buf
371
372 associate(u_vol => this%u_vol, v_vol => this%v_vol, &
373 w_vol => this%w_vol, p_vol => this%p_vol)
374
375 n = c_xh%dof%size()
376
377 ! If either dt or the backwards difference coefficient change,
378 ! then recompute base flow solution corresponding to unit forcing:
379
380 ifcomp = 0.0_rp
381
382 if ((.not. abscmp(dt, this%dtlag)) .or. &
383 (.not. abscmp(ext_bdf%diffusion_coeffs%x(1), this%bdlag))) then
384 ifcomp = 1.0_rp
385 end if
386
387 this%dtlag = dt
388 this%bdlag = ext_bdf%diffusion_coeffs%x(1)
389
390 call mpi_allreduce(mpi_in_place, ifcomp, 1, &
391 mpi_real_precision, mpi_sum, neko_comm, ierr)
392
393 if (ifcomp .gt. 0d0) then
394 call this%compute(u_res, v_res, w_res, p_res, &
395 ext_bdf, gs_xh, c_xh, rho, mu, ext_bdf%diffusion_coeffs%x(1), dt, &
396 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
397 ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
398 vel_max_iter)
399 end if
400
401 if (neko_bcknd_device .eq. 1) then
402 if (this%flow_dir .eq. 1) then
403 current_flow = &
404 device_glsc2(u%x_d, c_xh%B_d, n) / this%domain_length ! for X
405 else if (this%flow_dir .eq. 2) then
406 current_flow = &
407 device_glsc2(v%x_d, c_xh%B_d, n) / this%domain_length ! for Y
408 else if (this%flow_dir .eq. 3) then
409 current_flow = &
410 device_glsc2(w%x_d, c_xh%B_d, n) / this%domain_length ! for Z
411 end if
412 else
413 if (this%flow_dir .eq. 1) then
414 current_flow = glsc2(u%x, c_xh%B, n) / this%domain_length ! for X
415 else if (this%flow_dir .eq. 2) then
416 current_flow = glsc2(v%x, c_xh%B, n) / this%domain_length ! for Y
417 else if (this%flow_dir .eq. 3) then
418 current_flow = glsc2(w%x, c_xh%B, n) / this%domain_length ! for Z
419 end if
420 end if
421
422 if (this%avflow) then
423 xsec = c_xh%volume / this%domain_length
424 flow_rate = this%flow_rate*xsec
425 else
426 flow_rate = this%flow_rate
427 end if
428
429 delta_flow = flow_rate - current_flow
430 scale = delta_flow / this%base_flow
431
432 if (this%log .and. pe_rank .eq. 0) then
433 if (this%flow_dir .eq. 1) then
434 flow_dir_label = ' x'
435 else if (this%flow_dir .eq. 2) then
436 flow_dir_label = ' y'
437 else if (this%flow_dir .eq. 3) then
438 flow_dir_label = ' z'
439 end if
440 write(step_str, '(I12)') time%tstep
441 step_str = adjustl(step_str)
442 write(log_buf, '(A,A3,A5,1X,5A18)') &
443 'Flow rate ', ' | ', 'Dir.:', 'Time:', 'Scale:', 'Rate:', &
444 'Current:', 'Base:'
445 call neko_log%message(log_buf)
446 write(log_buf, '(A12,A3,A5,1X,5E18.9)') step_str, ' | ', &
447 flow_dir_label, time%t, scale, flow_rate, current_flow, &
448 this%base_flow
449 call neko_log%message(log_buf)
450 end if
451
452 if (neko_bcknd_device .eq. 1) then
453 call device_add2s2(u%x_d, u_vol%x_d, scale, n)
454 call device_add2s2(v%x_d, v_vol%x_d, scale, n)
455 call device_add2s2(w%x_d, w_vol%x_d, scale, n)
456 call device_add2s2(p%x_d, p_vol%x_d, scale, n)
457 else
458 do concurrent(i = 1: n)
459 u%x(i,1,1,1) = u%x(i,1,1,1) + scale * u_vol%x(i,1,1,1)
460 v%x(i,1,1,1) = v%x(i,1,1,1) + scale * v_vol%x(i,1,1,1)
461 w%x(i,1,1,1) = w%x(i,1,1,1) + scale * w_vol%x(i,1,1,1)
462 p%x(i,1,1,1) = p%x(i,1,1,1) + scale * p_vol%x(i,1,1,1)
463 end do
464 end if
465 end associate
466
467 end subroutine fluid_vol_flow
468
469end module fluid_volflow
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_opchsign(a1_d, a2_d, a3_d, gdim, n)
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
subroutine fluid_vol_flow(this, u, v, w, p, u_res, v_res, w_res, p_res, c_xh, gs_xh, ext_bdf, rho, mu, dt, time, bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, vel_max_iter)
Adjust flow volume.
subroutine fluid_vol_flow_init(this, dm_xh, params)
subroutine fluid_vol_flow_compute(this, u_res, v_res, w_res, p_res, ext_bdf, gs_xh, c_xh, rho, mu, bd, dt, bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, vel_max_iter)
Compute flow adjustment.
subroutine fluid_vol_flow_free(this)
Gather-scatter.
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition mathops.f90:65
subroutine, public opchsign(a1, a2, a3, gdim, n)
for and .
Definition mathops.f90:76
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
Krylov preconditioner.
Definition precon.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.
Compound scheme for the advection and diffusion operators in a transport equation.
Module with things related to the simulation time.
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:61
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
A struct that contains all info about the time, expand as needed.