Neko 1.99.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", &
217 ignore_existing = .true.)
218
219 call this%abx2%init(dm_xh, trim(this%name) // "_abx2")
220 call neko_registry%add_field(dm_xh, trim(this%name) // "_abx2", &
221 ignore_existing = .true.)
222
223 call this%advs%init(dm_xh, "advs")
224
225 call this%ds%init(dm_xh, 'ds')
226
227 end associate
228
229 ! Set up boundary conditions
230 call this%setup_bcs_(user)
231
232 ! Initialize dirichlet bcs for scalar residual
233 call this%bc_res%init(this%c_Xh, params)
234 do i = 1, this%bcs%size()
235 if (this%bcs%strong(i)) then
236 bc_i => this%bcs%get(i)
237 call this%bc_res%mark_facets(bc_i%marked_facet)
238 end if
239 end do
240
241! call this%bc_res%mark_zones_from_list('d_s', this%bc_labels)
242 call this%bc_res%finalize()
243
244 call this%bclst_ds%init()
245 call this%bclst_ds%append(this%bc_res)
246
247
248 ! Initialize projection space
249 call this%proj_s%init(this%dm_Xh%size(), this%projection_dim, &
250 this%projection_activ_step)
251
252 ! Determine the time-interpolation scheme
253 call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
254 ! Point to case checkpoint
255 this%chkp => chkp
256 ! Initialize advection factory
257 call json_get_or_default(params, 'advection', advection, .true.)
258
259 call advection_factory(this%adv, numerics_params, this%c_Xh, &
260 ulag, vlag, wlag, this%chkp%dtlag, &
261 this%chkp%tlag, time_scheme, .not. advection, &
262 this%slag)
263 end subroutine scalar_pnpn_init
264
265 ! Restarts the scalar from a checkpoint
266 subroutine scalar_pnpn_restart(this, chkp)
267 class(scalar_pnpn_t), target, intent(inout) :: this
268 type(chkp_t), intent(inout) :: chkp
269 real(kind=rp) :: dtlag(10), tlag(10)
270 integer :: n
271 type(field_t), pointer :: temp_field
272 dtlag = chkp%dtlag
273 tlag = chkp%tlag
274
275 n = this%s%dof%size()
276
277 ! Lag fields are restored through the checkpoint's fsp mechanism
278
279 call col2(this%s%x, this%c_Xh%mult, n)
280 call col2(this%slag%lf(1)%x, this%c_Xh%mult, n)
281 call col2(this%slag%lf(2)%x, this%c_Xh%mult, n)
282 if (neko_bcknd_device .eq. 1) then
283 call device_memcpy(this%s%x, this%s%x_d, &
284 n, host_to_device, sync = .false.)
285 call device_memcpy(this%slag%lf(1)%x, this%slag%lf(1)%x_d, &
286 n, host_to_device, sync = .false.)
287 call device_memcpy(this%slag%lf(2)%x, this%slag%lf(2)%x_d, &
288 n, host_to_device, sync = .false.)
289 call device_memcpy(this%abx1%x, this%abx1%x_d, &
290 n, host_to_device, sync = .false.)
291 call device_memcpy(this%abx2%x, this%abx2%x_d, &
292 n, host_to_device, sync = .false.)
293 call device_memcpy(this%advs%x, this%advs%x_d, &
294 n, host_to_device, sync = .false.)
295 end if
296
297 call this%gs_Xh%op(this%s, gs_op_add)
298 call this%gs_Xh%op(this%slag%lf(1), gs_op_add)
299 call this%gs_Xh%op(this%slag%lf(2), gs_op_add)
300
301 end subroutine scalar_pnpn_restart
302
303 subroutine scalar_pnpn_free(this)
304 class(scalar_pnpn_t), intent(inout) :: this
305
306 !Deallocate scalar field
307 call this%scheme_free()
308
309 call this%bc_res%free()
310 call this%bclst_ds%free()
311 call this%proj_s%free()
312
313 call this%s_res%free()
314
315 call this%ds%free()
316
317 call this%abx1%free()
318 call this%abx2%free()
319
320 call this%advs%free()
321
322 if (allocated(this%adv)) then
323 call this%adv%free()
324 deallocate(this%adv)
325 end if
326
327 if (allocated(this%Ax)) then
328 deallocate(this%Ax)
329 end if
330
331 if (allocated(this%res)) then
332 deallocate(this%res)
333 end if
334
335 if (allocated(this%makeext)) then
336 deallocate(this%makeext)
337 end if
338
339 if (allocated(this%makebdf)) then
340 deallocate(this%makebdf)
341 end if
342
343 if (allocated(this%makeoifs)) then
344 deallocate(this%makeoifs)
345 end if
346
347 end subroutine scalar_pnpn_free
348
349 subroutine scalar_pnpn_step(this, time, ext_bdf, dt_controller, &
350 ksp_results)
351 class(scalar_pnpn_t), intent(inout) :: this
352 type(time_state_t), intent(in) :: time
353 type(time_scheme_controller_t), intent(in) :: ext_bdf
354 type(time_step_controller_t), intent(in) :: dt_controller
355 type(ksp_monitor_t), intent(inout) :: ksp_results
356 ! Number of degrees of freedom
357 integer :: n
358
359 if (this%freeze) return
360
361 n = this%dm_Xh%size()
362
363 call profiler_start_region(trim(this%name), 2)
364 associate(u => this%u, v => this%v, w => this%w, s => this%s, &
365 cp => this%cp, rho => this%rho, lambda_tot => this%lambda_tot, &
366 ds => this%ds, &
367 s_res => this%s_res, &
368 ax => this%Ax, f_xh => this%f_Xh, xh => this%Xh, &
369 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
370 slag => this%slag, oifs => this%oifs, &
371 projection_dim => this%projection_dim, &
372 msh => this%msh, res => this%res, makeoifs => this%makeoifs, &
373 makeext => this%makeext, makebdf => this%makebdf, &
374 t => time%t, tstep => time%tstep, dt => time%dt)
375
376 ! Logs extra information the log level is NEKO_LOG_DEBUG or above.
377 call print_debug(this)
378 ! Compute the source terms
379 call this%source_term%compute(time)
380
381 ! Apply weak boundary conditions, that contribute to the source terms.
382 call this%bcs%apply_scalar(this%f_Xh%x, dm_xh%size(), time, .false.)
383
384 if (oifs) then
385 ! Add the advection operators to the right-hans-side.
386 call this%adv%compute_scalar(u, v, w, s, this%advs, &
387 xh, this%c_Xh, dm_xh%size())
388
389 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
390 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
391
392 call makeoifs%compute_scalar(this%advs%x, f_xh%x, rho%x(1,1,1,1), dt,&
393 n)
394 else
395 ! Add the advection operators to the right-hans-side.
396 call this%adv%compute_scalar(u, v, w, s, f_xh, &
397 xh, this%c_Xh, dm_xh%size())
398
399 ! At this point the RHS contains the sum of the advection operator,
400 ! Neumann boundary sources and additional source terms, evaluated using
401 ! the scalar field from the previous time-step. Now, this value is used in
402 ! the explicit time scheme to advance these terms in time.
403 call makeext%compute_scalar(this%abx1, this%abx2, f_xh%x, &
404 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
405
406 ! Add the RHS contributions coming from the BDF scheme.
407 call makebdf%compute_scalar(slag, f_xh%x, s, c_xh%B, rho%x(1,1,1,1), &
408 dt, ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
409 end if
410
411 call slag%update()
412
414 call this%apply_strong_bcs(time)
415
416 ! Update material properties if necessary
417 call this%update_material_properties(time)
418
419 ! Compute scalar residual.
420 call profiler_start_region(trim(this%name) // '_residual', 20)
421 call res%compute(ax, s, s_res, f_xh, c_xh, msh, xh, lambda_tot, &
422 rho%x(1,1,1,1)*cp%x(1,1,1,1), ext_bdf%diffusion_coeffs(1), dt, &
423 dm_xh%size())
424
425 call gs_xh%op(s_res, gs_op_add)
426
427 ! Apply a 0-valued Dirichlet boundary conditions on the ds.
428 call this%bclst_ds%apply_scalar(s_res%x, dm_xh%size())
429
430 call profiler_end_region(trim(this%name) // '_residual', 20)
431
432 call this%proj_s%pre_solving(s_res%x, tstep, c_xh, n, dt_controller)
433
434 call this%pc%update()
435 call profiler_start_region(trim(this%name) // '_solve', 21)
436 ksp_results = this%ksp%solve(ax, ds, s_res%x, n, &
437 c_xh, this%bclst_ds, gs_xh)
438 ksp_results%name = trim(this%name)
439 call profiler_end_region(trim(this%name) // '_solve', 21)
440
441 call this%proj_s%post_solving(ds%x, ax, c_xh, this%bclst_ds, gs_xh, &
442 n, tstep, dt_controller)
443
444 ! Update the solution
445 if (neko_bcknd_device .eq. 1) then
446 call device_add2s2(s%x_d, ds%x_d, 1.0_rp, n)
447 else
448 call add2s2(s%x, ds%x, 1.0_rp, n)
449 end if
450
451 end associate
452 call profiler_end_region(trim(this%name), 2)
453 end subroutine scalar_pnpn_step
454
455 subroutine print_debug(this)
456 class(scalar_pnpn_t), intent(inout) :: this
457 character(len=LOG_SIZE) :: log_buf
458 integer :: n
459
460 n = this%dm_Xh%size()
461
462 write(log_buf, '(A, A, E15.7, A, E15.7, A, E15.7)') 'Scalar debug', &
463 ' l2norm s', glsc2(this%s%x, this%s%x, n), &
464 ' slag1', glsc2(this%slag%lf(1)%x, this%slag%lf(1)%x, n), &
465 ' slag2', glsc2(this%slag%lf(2)%x, this%slag%lf(2)%x, n)
466 call neko_log%message(log_buf, lvl = neko_log_debug)
467 write(log_buf, '(A, A, E15.7, A, E15.7)') 'Scalar debug2', &
468 ' l2norm abx1', glsc2(this%abx1%x, this%abx1%x, n), &
469 ' abx2', glsc2(this%abx2%x, this%abx2%x, n)
470 call neko_log%message(log_buf, lvl = neko_log_debug)
471 end subroutine print_debug
472
475 subroutine scalar_pnpn_setup_bcs_(this, user)
476 class(scalar_pnpn_t), intent(inout) :: this
477 type(user_t), target, intent(in) :: user
478 integer :: i, j, n_bcs, zone_size, global_zone_size, ierr
479 type(json_core) :: core
480 type(json_value), pointer :: bc_object
481 type(json_file) :: bc_subdict
482 class(bc_t), pointer :: bc_i
483 logical :: found
484 ! Monitor which boundary zones have been marked
485 logical, allocatable :: marked_zones(:)
486 integer, allocatable :: zone_indices(:)
487
488 if (this%params%valid_path('boundary_conditions')) then
489 call this%params%info('boundary_conditions', &
490 n_children = n_bcs)
491 call this%params%get_core(core)
492 call this%params%get('boundary_conditions', bc_object, found)
493
494 call this%bcs%init(n_bcs)
495
496 allocate(marked_zones(size(this%msh%labeled_zones)))
497 marked_zones = .false.
498
499 do i = 1, n_bcs
500 ! Create a new json containing just the subdict for this bc
501 call json_extract_item(core, bc_object, i, bc_subdict)
502
503 ! Check that we are not trying to assing a bc to zone, for which one
504 ! has already been assigned and that the zone has more than 0 size
505 ! in the mesh.
506 call json_get(bc_subdict, "zone_indices", zone_indices)
507
508 do j = 1, size(zone_indices)
509 zone_size = this%msh%labeled_zones(zone_indices(j))%size
510 call mpi_allreduce(zone_size, global_zone_size, 1, &
511 mpi_integer, mpi_max, neko_comm, ierr)
512
513 if (global_zone_size .eq. 0) then
514 write(error_unit, '(A, A, I0, A, A, I0, A)') &
515 "*** ERROR ***: ", "Zone index ", zone_indices(j), &
516 " is invalid as this zone has 0 size, meaning it ", &
517 "does not exist in the mesh. Check scalar boundary ", &
518 "condition ", i, "."
519 error stop
520 end if
521
522 if (marked_zones(zone_indices(j))) then
523 write(error_unit, '(A, A, I0, A, A, A, A)') "*** ERROR ***: ", &
524 "Zone with index ", zone_indices(j), &
525 " has already been assigned a boundary condition. ", &
526 "Please check your boundary_conditions entry for the ", &
527 "scalar and make sure that each zone index appears only ",&
528 "in a single boundary condition."
529 error stop
530 else
531 marked_zones(zone_indices(j)) = .true.
532 end if
533 end do
534
535 bc_i => null()
536
537 call bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
538 call this%bcs%append(bc_i)
539 end do
540
541 ! Make sure all labeled zones with non-zero size have been marked
542 do i = 1, size(this%msh%labeled_zones)
543 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
544 (.not. marked_zones(i))) then
545 write(error_unit, '(A, A, I0)') "*** ERROR ***: ", &
546 "No scalar boundary condition assigned to zone ", i
547 error stop
548 end if
549 end do
550 else
551 ! Check that there are no labeled zones, i.e. all are periodic.
552 do i = 1, size(this%msh%labeled_zones)
553 if (this%msh%labeled_zones(i)%size .gt. 0) then
554 write(error_unit, '(A, A, A)') "*** ERROR ***: ", &
555 "No boundary_conditions entry in the case file for scalar ", &
556 this%s%name
557 error stop
558 end if
559 end do
560
561 ! For a pure periodic case, we still need to initilise the bc lists
562 ! to a zero size to avoid issues with apply() in step()
563 call this%bcs%init()
564
565 end if
566 end subroutine scalar_pnpn_setup_bcs_
567
570 subroutine scalar_scheme_apply_strong_bcs(this, time)
571 class(scalar_pnpn_t), intent(inout) :: this
572 type(time_state_t), intent(in) :: time
573
574 integer :: i
575 class(bc_t), pointer :: bc_i
576 bc_i => null()
577
578 ! First apply call, sets the Dirichlet value, let's call it d.
579 call this%bcs%apply(this%s, time = time, strong = .true.)
580 ! If we now have local nodes sharing the same global node, and with
581 ! some nodes not masked as Dirichlet, the node which *is* masked
582 ! will have the value d, and the the other ones just some value u.
583 ! Take a nodewise minimum between the local nodes.
584 ! Now, all local nodes store m = min(d, u)
585 call this%gs_Xh%op(this%s, gs_op_min, glb_cmd_event)
586 call device_event_sync(glb_cmd_event)
587
588 ! Second apply call, so Dirichlet nodes again store d, the rest still store
589 ! m, where m < d by construction.
590 call this%bcs%apply(this%s, time = time, strong = .true.)
591 ! Now apply a max, which guarantees that d wins and gets stored in all the
592 ! local nodes.
593 call this%gs_Xh%op(this%s, gs_op_max, glb_cmd_event)
594 call device_event_sync(glb_cmd_event)
595
596 ! Reset updated flags
597 do i = 1, this%bcs%size()
598 bc_i => this%bcs%get(i)
599 bc_i%updated = .false.
600 end do
601 nullify(bc_i)
602
603 end subroutine scalar_scheme_apply_strong_bcs
604
605
606end 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:87
type(log_t), public neko_log
Global log stream.
Definition log.f90:77
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:149
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:62
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:56
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...