Neko 1.0.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_pnpn.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2024, 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!
34
36 use num_types, only: rp
37 use, intrinsic :: iso_fortran_env, only: error_unit
39 rhs_maker_ext_fctry, rhs_maker_bdf_fctry, rhs_maker_oifs_fctry
41 use checkpoint, only : chkp_t
42 use field, only : field_t
43 use bc_list, only : bc_list_t
44 use mesh, only : mesh_t
45 use coefs, only : coef_t
48 use gather_scatter, only : gs_t, gs_op_add, gs_op_min, gs_op_max
49 use scalar_residual, only : scalar_residual_t, scalar_residual_factory
50 use ax_product, only : ax_t, ax_helm_factory
52 use registry, only: neko_registry
54 use krylov, only : ksp_monitor_t
57 use projection, only : projection_t
58 use math, only : glsc2, col2, add2s2
60 use advection, only : advection_t, advection_factory
63 use json_module, only : json_file, json_core, json_value
64 use user_intf, only : user_t
68 use time_state, only : time_state_t
69 use bc, only : bc_t
70 use comm, only : neko_comm
71 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_max
72 implicit none
73 private
74
75
76 type, public, extends(scalar_scheme_t) :: scalar_pnpn_t
77
79 type(field_t) :: s_res
80
82 type(field_t) :: ds
83
85 class(ax_t), allocatable :: ax
86
88 type(projection_t) :: proj_s
89
93 type(zero_dirichlet_t) :: bc_res
94
99 type(bc_list_t) :: bclst_ds
100
102 class(advection_t), allocatable :: adv
103
104 ! Time interpolation scheme
105 logical :: oifs
106
107 ! Advection terms for the oifs method
108 type(field_t) :: advs
109
111 class(scalar_residual_t), allocatable :: res
112
114 class(rhs_maker_ext_t), allocatable :: makeext
115
117 class(rhs_maker_bdf_t), allocatable :: makebdf
118
120 class(rhs_maker_oifs_t), allocatable :: makeoifs
121
123 type(field_t) :: abx1, abx2
124
125 contains
127 procedure, pass(this) :: init => scalar_pnpn_init
129 procedure, pass(this) :: restart => scalar_pnpn_restart
131 procedure, pass(this) :: free => scalar_pnpn_free
133 procedure, pass(this) :: step => scalar_pnpn_step
135 procedure, pass(this) :: apply_strong_bcs => scalar_scheme_apply_strong_bcs
137 procedure, pass(this) :: setup_bcs_ => scalar_pnpn_setup_bcs_
139 end type scalar_pnpn_t
140
141 interface
142
148 module subroutine bc_factory(object, scheme, json, coef, user)
149 class(bc_t), pointer, intent(inout) :: object
150 type(scalar_pnpn_t), intent(in) :: scheme
151 type(json_file), intent(inout) :: json
152 type(coef_t), intent(in) :: coef
153 type(user_t), intent(in) :: user
154 end subroutine bc_factory
155 end interface
156
157contains
158
171 subroutine scalar_pnpn_init(this, msh, coef, gs, params, numerics_params, &
172 user, chkp, ulag, vlag, wlag, time_scheme, rho)
173 class(scalar_pnpn_t), target, intent(inout) :: this
174 type(mesh_t), target, intent(in) :: msh
175 type(coef_t), target, intent(in) :: coef
176 type(gs_t), target, intent(inout) :: gs
177 type(json_file), target, intent(inout) :: params
178 type(json_file), target, intent(inout) :: numerics_params
179 type(user_t), target, intent(in) :: user
180 type(chkp_t), target, intent(inout) :: chkp
181 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
182 type(time_scheme_controller_t), target, intent(in) :: time_scheme
183 type(field_t), target, intent(in) :: rho
184 integer :: i
185 class(bc_t), pointer :: bc_i
186 character(len=15), parameter :: scheme = 'Modular (Pn/Pn)'
187 logical :: advection
188
189 call this%free()
190
191 ! Initiliaze base type.
192 call this%scheme_init(msh, coef, gs, params, scheme, user, rho)
193
194 ! Setup backend dependent Ax routines
195 call ax_helm_factory(this%ax, full_formulation = .false.)
196
197 ! Setup backend dependent scalar residual routines
198 call scalar_residual_factory(this%res)
199
200 ! Setup backend dependent summation of extrapolation scheme
201 call rhs_maker_ext_fctry(this%makeext)
202
203 ! Setup backend dependent contributions to F from lagged BD terms
204 call rhs_maker_bdf_fctry(this%makebdf)
205
206 ! Setup backend dependent contributions of the OIFS scheme
207 call rhs_maker_oifs_fctry(this%makeoifs)
208
209 ! Initialize variables specific to this plan
210 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
211 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
212
213 call this%s_res%init(dm_xh, "s_res")
214
215 call this%abx1%init(dm_xh, trim(this%name)//"_abx1")
216 call neko_registry%add_field(dm_xh, trim(this%name)//"_abx1", ignore_existing = .true.)
217
218 call this%abx2%init(dm_xh, trim(this%name)//"_abx2")
219 call neko_registry%add_field(dm_xh, trim(this%name)//"_abx2", ignore_existing = .true.)
220
221 call this%advs%init(dm_xh, "advs")
222
223 call this%ds%init(dm_xh, 'ds')
224
225 end associate
226
227 ! Set up boundary conditions
228 call this%setup_bcs_(user)
229
230 ! Initialize dirichlet bcs for scalar residual
231 call this%bc_res%init(this%c_Xh, params)
232 do i = 1, this%bcs%size()
233 if (this%bcs%strong(i)) then
234 bc_i => this%bcs%get(i)
235 call this%bc_res%mark_facets(bc_i%marked_facet)
236 end if
237 end do
238
239! call this%bc_res%mark_zones_from_list('d_s', this%bc_labels)
240 call this%bc_res%finalize()
241
242 call this%bclst_ds%init()
243 call this%bclst_ds%append(this%bc_res)
244
245
246 ! Initialize projection space
247 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
248 this%projection_activ_step)
249
250 ! Determine the time-interpolation scheme
251 call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
252 ! Point to case checkpoint
253 this%chkp => chkp
254 ! Initialize advection factory
255 call json_get_or_default(params, 'advection', advection, .true.)
256
257 call advection_factory(this%adv, numerics_params, this%c_Xh, &
258 ulag, vlag, wlag, this%chkp%dtlag, &
259 this%chkp%tlag, time_scheme, .not. advection, &
260 this%slag)
261 end subroutine scalar_pnpn_init
262
263 ! Restarts the scalar from a checkpoint
264 subroutine scalar_pnpn_restart(this, chkp)
265 class(scalar_pnpn_t), target, intent(inout) :: this
266 type(chkp_t), intent(inout) :: chkp
267 real(kind=rp) :: dtlag(10), tlag(10)
268 integer :: n
269 type(field_t), pointer :: temp_field
270 dtlag = chkp%dtlag
271 tlag = chkp%tlag
272
273 n = this%s%dof%size()
274
275 ! Lag fields are restored through the checkpoint's fsp mechanism
276
277 call col2(this%s%x, this%c_Xh%mult, n)
278 call col2(this%slag%lf(1)%x, this%c_Xh%mult, n)
279 call col2(this%slag%lf(2)%x, this%c_Xh%mult, n)
280 if (neko_bcknd_device .eq. 1) then
281 call device_memcpy(this%s%x, this%s%x_d, &
282 n, host_to_device, sync = .false.)
283 call device_memcpy(this%slag%lf(1)%x, this%slag%lf(1)%x_d, &
284 n, host_to_device, sync = .false.)
285 call device_memcpy(this%slag%lf(2)%x, this%slag%lf(2)%x_d, &
286 n, host_to_device, sync = .false.)
287 call device_memcpy(this%abx1%x, this%abx1%x_d, &
288 n, host_to_device, sync = .false.)
289 call device_memcpy(this%abx2%x, this%abx2%x_d, &
290 n, host_to_device, sync = .false.)
291 call device_memcpy(this%advs%x, this%advs%x_d, &
292 n, host_to_device, sync = .false.)
293 end if
294
295 call this%gs_Xh%op(this%s, gs_op_add)
296 call this%gs_Xh%op(this%slag%lf(1), gs_op_add)
297 call this%gs_Xh%op(this%slag%lf(2), gs_op_add)
298
299 end subroutine scalar_pnpn_restart
300
301 subroutine scalar_pnpn_free(this)
302 class(scalar_pnpn_t), intent(inout) :: this
303
304 !Deallocate scalar field
305 call this%scheme_free()
306
307 call this%bc_res%free()
308 call this%bclst_ds%free()
309 call this%proj_s%free()
310
311 call this%s_res%free()
312
313 call this%ds%free()
314
315 call this%abx1%free()
316 call this%abx2%free()
317
318 call this%advs%free()
319
320 if (allocated(this%adv)) then
321 call this%adv%free()
322 deallocate(this%adv)
323 end if
324
325 if (allocated(this%Ax)) then
326 deallocate(this%Ax)
327 end if
328
329 if (allocated(this%res)) then
330 deallocate(this%res)
331 end if
332
333 if (allocated(this%makeext)) then
334 deallocate(this%makeext)
335 end if
336
337 if (allocated(this%makebdf)) then
338 deallocate(this%makebdf)
339 end if
340
341 if (allocated(this%makeoifs)) then
342 deallocate(this%makeoifs)
343 end if
344
345 end subroutine scalar_pnpn_free
346
347 subroutine scalar_pnpn_step(this, time, ext_bdf, dt_controller, &
348 ksp_results)
349 class(scalar_pnpn_t), intent(inout) :: this
350 type(time_state_t), intent(in) :: time
351 type(time_scheme_controller_t), intent(in) :: ext_bdf
352 type(time_step_controller_t), intent(in) :: dt_controller
353 type(ksp_monitor_t), intent(inout) :: ksp_results
354 ! Number of degrees of freedom
355 integer :: n
356
357 if (this%freeze) return
358
359 n = this%dm_Xh%size()
360
361 call profiler_start_region(trim(this%name), 2)
362 associate(u => this%u, v => this%v, w => this%w, s => this%s, &
363 cp => this%cp, rho => this%rho, lambda_tot => this%lambda_tot, &
364 ds => this%ds, &
365 s_res => this%s_res, &
366 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
367 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
368 slag => this%slag, oifs => this%oifs, &
369 projection_dim => this%projection_dim, &
370 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
371 makeext => this%makeext, makebdf => this%makebdf, &
372 t => time%t, tstep => time%tstep, dt => time%dt)
373
374 ! Logs extra information the log level is NEKO_LOG_DEBUG or above.
375 call print_debug(this)
376 ! Compute the source terms
377 call this%source_term%compute(time)
378
379 ! Apply weak boundary conditions, that contribute to the source terms.
380 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), time, .false.)
381
382 if (oifs) then
383 ! Add the advection operators to the right-hans-side.
384 call this%adv%compute_scalar(u, v, w, s, this%advs, &
385 xh, this%c_Xh, dm_xh%size())
386
387 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
388 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
389
390 call makeoifs%compute_scalar(this%advs%x, f_xh%x, rho%x(1,1,1,1), dt,&
391 n)
392 else
393 ! Add the advection operators to the right-hans-side.
394 call this%adv%compute_scalar(u, v, w, s, f_xh, &
395 xh, this%c_Xh, dm_xh%size())
396
397 ! At this point the RHS contains the sum of the advection operator,
398 ! Neumann boundary sources and additional source terms, evaluated using
399 ! the scalar field from the previous time-step. Now, this value is used in
400 ! the explicit time scheme to advance these terms in time.
401 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
402 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
403
404 ! Add the RHS contributions coming from the BDF scheme.
405 call makebdf%compute_scalar(slag, f_xh%x, s, c_xh%B, rho%x(1,1,1,1), &
406 dt, ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
407 end if
408
409 call slag%update()
410
412 call this%apply_strong_bcs(time)
413
414 ! Update material properties if necessary
415 call this%update_material_properties(time)
416
417 ! Compute scalar residual.
418 call profiler_start_region(trim(this%name) // '_residual', 20)
419 call res%compute(ax, s, s_res, f_xh, c_xh, msh, xh, lambda_tot, &
420 rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs(1), dt, &
421 dm_xh%size())
422
423 call gs_xh%op(s_res, gs_op_add)
424
425 ! Apply a 0-valued Dirichlet boundary conditions on the ds.
426 call this%bclst_ds%apply_scalar(s_res%x, dm_xh%size())
427
428 call profiler_end_region(trim(this%name) // '_residual', 20)
429
430 call this%proj_s%pre_solving(s_res%x, tstep, c_xh, n, dt_controller)
431
432 call this%pc%update()
433 call profiler_start_region(trim(this%name) // '_solve', 21)
434 ksp_results = this%ksp%solve(ax, ds, s_res%x, n, &
435 c_xh, this%bclst_ds, gs_xh)
436 ksp_results%name = trim(this%name)
437 call profiler_end_region(trim(this%name) // '_solve', 21)
438
439 call this%proj_s%post_solving(ds%x, ax, c_xh, this%bclst_ds, gs_xh, &
440 n, tstep, dt_controller)
441
442 ! Update the solution
443 if (neko_bcknd_device .eq. 1) then
444 call device_add2s2(s%x_d, ds%x_d, 1.0_rp, n)
445 else
446 call add2s2(s%x, ds%x, 1.0_rp, n)
447 end if
448
449 end associate
450 call profiler_end_region(trim(this%name), 2)
451 end subroutine scalar_pnpn_step
452
453 subroutine print_debug(this)
454 class(scalar_pnpn_t), intent(inout) :: this
455 character(len=LOG_SIZE) :: log_buf
456 integer :: n
457
458 n = this%dm_Xh%size()
459
460 write(log_buf, '(A, A, E15.7, A, E15.7, A, E15.7)') 'Scalar debug', &
461 ' l2norm s', glsc2(this%s%x, this%s%x, n), &
462 ' slag1', glsc2(this%slag%lf(1)%x, this%slag%lf(1)%x, n), &
463 ' slag2', glsc2(this%slag%lf(2)%x, this%slag%lf(2)%x, n)
464 call neko_log%message(log_buf, lvl = neko_log_debug)
465 write(log_buf, '(A, A, E15.7, A, E15.7)') 'Scalar debug2', &
466 ' l2norm abx1', glsc2(this%abx1%x, this%abx1%x, n), &
467 ' abx2', glsc2(this%abx2%x, this%abx2%x, n)
468 call neko_log%message(log_buf, lvl = neko_log_debug)
469 end subroutine print_debug
470
473 subroutine scalar_pnpn_setup_bcs_(this, user)
474 class(scalar_pnpn_t), intent(inout) :: this
475 type(user_t), target, intent(in) :: user
476 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
477 type(json_core) :: core
478 type(json_value), pointer :: bc_object
479 type(json_file) :: bc_subdict
480 class(bc_t), pointer :: bc_i
481 logical :: found
482 ! Monitor which boundary zones have been marked
483 logical, allocatable :: marked_zones(:)
484 integer, allocatable :: zone_indices(:)
485
486 if (this%params%valid_path('boundary_conditions')) then
487 call this%params%info('boundary_conditions', &
488 n_children = n_bcs)
489 call this%params%get_core(core)
490 call this%params%get('boundary_conditions', bc_object, found)
491
492 call this%bcs%init(n_bcs)
493
494 allocate(marked_zones(size(this%msh%labeled_zones)))
495 marked_zones = .false.
496
497 do i = 1, n_bcs
498 ! Create a new json containing just the subdict for this bc
499 call json_extract_item(core, bc_object, i, bc_subdict)
500
501 ! Check that we are not trying to assing a bc to zone, for which one
502 ! has already been assigned and that the zone has more than 0 size
503 ! in the mesh.
504 call json_get(bc_subdict, "zone_indices", zone_indices)
505
506 do j = 1, size(zone_indices)
507 zone_size = this%msh%labeled_zones(zone_indices(j))%size
508 call mpi_allreduce(zone_size, global_zone_size, 1, &
509 mpi_integer, mpi_max, neko_comm, ierr)
510
511 if (global_zone_size .eq. 0) then
512 write(error_unit, '(A, A, I0, A, A, I0, A)') "*** ERROR ***: ",&
513 "Zone index ", zone_indices(j), &
514 " is invalid as this zone has 0 size, meaning it ", &
515 "does not exist in the mesh. Check scalar boundary condition ", &
516 i, "."
517 error stop
518 end if
519
520 if (marked_zones(zone_indices(j)) .eqv. .true.) then
521 write(error_unit, '(A, A, I0, A, A, A, A)') "*** ERROR ***: ", &
522 "Zone with index ", zone_indices(j), &
523 " has already been assigned a boundary condition. ", &
524 "Please check your boundary_conditions entry for the ", &
525 "scalar and make sure that each zone index appears only ",&
526 "in a single boundary condition."
527 error stop
528 else
529 marked_zones(zone_indices(j)) = .true.
530 end if
531 end do
532
533 bc_i => null()
534
535 call bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
536 call this%bcs%append(bc_i)
537 end do
538
539 ! Make sure all labeled zones with non-zero size have been marked
540 do i = 1, size(this%msh%labeled_zones)
541 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
542 (marked_zones(i) .eqv. .false.)) then
543 write(error_unit, '(A, A, I0)') "*** ERROR ***: ", &
544 "No scalar boundary condition assigned to zone ", i
545 error stop
546 end if
547 end do
548 else
549 ! Check that there are no labeled zones, i.e. all are periodic.
550 do i = 1, size(this%msh%labeled_zones)
551 if (this%msh%labeled_zones(i)%size .gt. 0) then
552 write(error_unit, '(A, A, A)') "*** ERROR ***: ", &
553 "No boundary_conditions entry in the case file for scalar ", &
554 this%s%name
555 error stop
556 end if
557 end do
558
559 ! For a pure periodic case, we still need to initilise the bc lists
560 ! to a zero size to avoid issues with apply() in step()
561 call this%bcs%init()
562
563 end if
564 end subroutine scalar_pnpn_setup_bcs_
565
568 subroutine scalar_scheme_apply_strong_bcs(this, time)
569 class(scalar_pnpn_t), intent(inout) :: this
570 type(time_state_t), intent(in) :: time
571
572 integer :: i
573 class(bc_t), pointer :: bc_i
574 bc_i => null()
575
576 ! First apply call, sets the Dirichlet value, let's call it d.
577 call this%bcs%apply(this%s, time = time, strong = .true.)
578 ! If we now have local nodes sharing the same global node, and with
579 ! some nodes not masked as Dirichlet, the node which *is* masked
580 ! will have the value d, and the the other ones just some value u.
581 ! Take a nodewise minimum between the local nodes.
582 ! Now, all local nodes store m = min(d, u)
583 call this%gs_Xh%op(this%s, gs_op_min, glb_cmd_event)
584 call device_event_sync(glb_cmd_event)
585
586 ! Second apply call, so Dirichlet nodes again store d, the rest still store
587 ! m, where m < d by construction.
588 call this%bcs%apply(this%s, time = time, strong = .true.)
589 ! Now apply a max, which guarantees that d wins and gets stored in all the
590 ! local nodes.
591 call this%gs_Xh%op(this%s, gs_op_max, glb_cmd_event)
592 call device_event_sync(glb_cmd_event)
593
594 ! Reset updated flags
595 do i = 1, this%bcs%size()
596 bc_i => this%bcs%get(i)
597 bc_i%updated = .false.
598 end do
599 nullify(bc_i)
600
601 end subroutine scalar_scheme_apply_strong_bcs
602
603
604end module scalar_pnpn
Copy data between host and device (or device and device)
Definition device.F90:71
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.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
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_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_event_sync(event)
Synchronize an event.
Definition device.F90:1314
integer, parameter, public host_to_device
Definition device.F90:47
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:62
Dirichlet condition applied in the facet normal direction.
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
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
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:84
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
Project x onto X, the space of old solutions and back again.
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:128
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition rhs_maker.f90:38
Contains the scalar_pnpn_t type.
subroutine scalar_pnpn_step(this, time, ext_bdf, dt_controller, ksp_results)
subroutine scalar_pnpn_setup_bcs_(this, user)
Initialize boundary conditions.
subroutine scalar_scheme_apply_strong_bcs(this, time)
Apply strong boundary conditions.
subroutine scalar_pnpn_restart(this, chkp)
subroutine scalar_pnpn_init(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Boundary condition factory. Both constructs and initializes the object.
subroutine scalar_pnpn_free(this)
Defines the residual for the scalar transport equation.
Contains the scalar_scheme_t type.
Compound scheme for the advection and diffusion operators in a transport equation.
Base class for time integration schemes.
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Defines a zero-valued Dirichlet boundary condition.
Base abstract type for computing the advection operator.
Definition advection.f90:46
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:61
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
Dirichlet condition in facet normal direction.
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Abstract type to add contributions to F from lagged BD terms.
Definition rhs_maker.f90:59
Abstract type to sum up contributions to kth order extrapolation scheme.
Definition rhs_maker.f90:52
Abstract type to add contributions of kth order OIFS scheme.
Definition rhs_maker.f90:66
Abstract type to compute scalar residual.
Base type for a scalar advection-diffusion solver.
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.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...