Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_pnpn.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, 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!
35 use, intrinsic :: iso_fortran_env, only: error_unit
36 use coefs, only : coef_t
37 use symmetry, only : symmetry_t
39 use logger, only: neko_log, log_size
40 use num_types, only : rp
41 use krylov, only : ksp_monitor_t
43 pnpn_prs_res_factory, pnpn_vel_res_factory, &
44 pnpn_prs_res_stress_factory, pnpn_vel_res_stress_factory
46 rhs_maker_oifs_t, rhs_maker_sumab_fctry, rhs_maker_bdf_fctry, &
47 rhs_maker_ext_fctry, rhs_maker_oifs_fctry
51 use fluid_aux, only : fluid_step_info
52 use projection, only : projection_t
56 use advection, only : advection_t, advection_factory
58 use json_module, only : json_file, json_core, json_value
60 use json_module, only : json_file
61 use ax_product, only : ax_t, ax_helm_factory
62 use field, only : field_t
63 use dirichlet, only : dirichlet_t
67 use non_normal, only : non_normal_t
68 use checkpoint, only : chkp_t
69 use mesh, only : mesh_t
70 use user_intf, only : user_t
72 use gs_ops, only : gs_op_add
74 use mathops, only : opadd2cm, opcolv
75 use bc_list, only: bc_list_t
79 use bc, only : bc_t
80 use file, only : file_t
81 use operators, only : ortho
82 use time_state, only : time_state_t
83 use comm, only : neko_comm
84 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_max, mpi_lor, &
85 mpi_integer, mpi_logical
86 implicit none
87 private
88
89
90
92
94 type(field_t) :: p_res, u_res, v_res, w_res
95
98 type(field_t) :: dp, du, dv, dw
99
100 !
101 ! Implicit operators, i.e. the left-hand-side of the Helmholz problem.
102 !
103
104 ! Coupled Helmholz operator for velocity
105 class(ax_t), allocatable :: ax_vel
106 ! Helmholz operator for pressure
107 class(ax_t), allocatable :: ax_prs
108
109 !
110 ! Projections for solver speed-up
111 !
112
114 type(projection_t) :: proj_prs
115 type(projection_vel_t) :: proj_vel
116
117 !
118 ! Special Karniadakis scheme boundary conditions in the pressure equation
119 !
120
122 type(facet_normal_t) :: bc_prs_surface
123
125 type(facet_normal_t) :: bc_sym_surface
126
127 !
128 ! Boundary conditions and lists for residuals and solution increments
129 !
130
132 type(zero_dirichlet_t) :: bc_vel_res
134 type(zero_dirichlet_t) :: bc_du
136 type(zero_dirichlet_t) :: bc_dv
138 type(zero_dirichlet_t) :: bc_dw
140 type(zero_dirichlet_t) :: bc_dp
141
143 type(bc_list_t) :: bclst_vel_res
144 type(bc_list_t) :: bclst_du
145 type(bc_list_t) :: bclst_dv
146 type(bc_list_t) :: bclst_dw
147 type(bc_list_t) :: bclst_dp
148
149
150 ! Checker for wether we have a strong pressure bc. If not, the pressure
151 ! is demeaned at every time step.
152 logical :: prs_dirichlet = .false.
153
154
155 ! The advection operator.
156 class(advection_t), allocatable :: adv
157
158 ! Time OIFS interpolation scheme for advection.
159 logical :: oifs
160
161 ! Time variables
162 type(field_t) :: abx1, aby1, abz1
163 type(field_t) :: abx2, aby2, abz2
164
165 ! Advection terms for the oifs method
166 type(field_t) :: advx, advy, advz
167
169 class(pnpn_prs_res_t), allocatable :: prs_res
170
172 class(pnpn_vel_res_t), allocatable :: vel_res
173
175 class(rhs_maker_sumab_t), allocatable :: sumab
176
178 class(rhs_maker_ext_t), allocatable :: makeabf
179
181 class(rhs_maker_bdf_t), allocatable :: makebdf
182
184 class(rhs_maker_oifs_t), allocatable :: makeoifs
185
187 type(fluid_volflow_t) :: vol_flow
188
190 logical :: full_stress_formulation = .false.
191
192 contains
194 procedure, pass(this) :: init => fluid_pnpn_init
196 procedure, pass(this) :: free => fluid_pnpn_free
198 procedure, pass(this) :: step => fluid_pnpn_step
200 procedure, pass(this) :: restart => fluid_pnpn_restart
202 procedure, pass(this) :: setup_bcs => fluid_pnpn_setup_bcs
204 procedure, pass(this) :: write_boundary_conditions => &
206 end type fluid_pnpn_t
207
208 interface
209
216 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
217 class(bc_t), pointer, intent(inout) :: object
218 type(fluid_pnpn_t), intent(in) :: scheme
219 type(json_file), intent(inout) :: json
220 type(coef_t), target, intent(in) :: coef
221 type(user_t), intent(in) :: user
222 end subroutine pressure_bc_factory
223 end interface
224
225 interface
226
233 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
234 class(bc_t), pointer, intent(inout) :: object
235 type(fluid_pnpn_t), intent(inout) :: scheme
236 type(json_file), intent(inout) :: json
237 type(coef_t), target, intent(in) :: coef
238 type(user_t), intent(in) :: user
239 end subroutine velocity_bc_factory
240 end interface
241
242contains
243
244 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
245 class(fluid_pnpn_t), target, intent(inout) :: this
246 type(mesh_t), target, intent(inout) :: msh
247 integer, intent(in) :: lx
248 type(json_file), target, intent(inout) :: params
249 type(user_t), target, intent(in) :: user
250 type(chkp_t), target, intent(inout) :: chkp
251 character(len=15), parameter :: scheme = 'Modular (Pn/Pn)'
252 integer :: i
253 class(bc_t), pointer :: bc_i, vel_bc
254 real(kind=rp) :: abs_tol
255 character(len=LOG_SIZE) :: log_buf
256 integer :: ierr, integer_val, solver_maxiter
257 character(len=:), allocatable :: solver_type, precon_type
258 logical :: monitor, found
259 logical :: advection
260 type(json_file) :: numerics_params, precon_params
261
262 call this%free()
263
264 ! Initialize base class
265 call this%init_base(msh, lx, params, scheme, user, .true.)
266
267 ! Add pressure field to the registry. For this scheme it is in the same
268 ! Xh as the velocity
269 call neko_field_registry%add_field(this%dm_Xh, 'p')
270 this%p => neko_field_registry%get_field('p')
271
272 !
273 ! Select governing equations via associated residual and Ax types
274 !
275
276 call json_get(params, 'case.numerics.time_order', integer_val)
277 allocate(this%ext_bdf)
278 call this%ext_bdf%init(integer_val)
279
280 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
281 this%full_stress_formulation, .false.)
282
283 if (this%full_stress_formulation .eqv. .true.) then
284 ! Setup backend dependent Ax routines
285 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
286
287 ! Setup backend dependent prs residual routines
288 call pnpn_prs_res_stress_factory(this%prs_res)
289
290 ! Setup backend dependent vel residual routines
291 call pnpn_vel_res_stress_factory(this%vel_res)
292 else
293 ! Setup backend dependent Ax routines
294 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
295
296 ! Setup backend dependent prs residual routines
297 call pnpn_prs_res_factory(this%prs_res)
298
299 ! Setup backend dependent vel residual routines
300 call pnpn_vel_res_factory(this%vel_res)
301 end if
302
303
304
305 if (params%valid_path('case.fluid.nut_field')) then
306 if (this%full_stress_formulation .eqv. .false.) then
307 call neko_error("You need to set full_stress_formulation to " // &
308 "true for the fluid to have a spatially varying " // &
309 "viscocity field.")
310 end if
311 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
312 else
313 this%nut_field_name = ""
314 end if
315
316 ! Setup Ax for the pressure
317 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
318
319
320 ! Setup backend dependent summation of AB/BDF
321 call rhs_maker_sumab_fctry(this%sumab)
322
323 ! Setup backend dependent summation of extrapolation scheme
324 call rhs_maker_ext_fctry(this%makeabf)
325
326 ! Setup backend depenent contributions to F from lagged BD terms
327 call rhs_maker_bdf_fctry(this%makebdf)
328
329 ! Setup backend dependent summations of the OIFS method
330 call rhs_maker_oifs_fctry(this%makeoifs)
331
332 ! Initialize variables specific to this plan
333 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
334 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
335
336 call this%p_res%init(dm_xh, "p_res")
337 call this%u_res%init(dm_xh, "u_res")
338 call this%v_res%init(dm_xh, "v_res")
339 call this%w_res%init(dm_xh, "w_res")
340 call this%abx1%init(dm_xh, "abx1")
341 call this%aby1%init(dm_xh, "aby1")
342 call this%abz1%init(dm_xh, "abz1")
343 call this%abx2%init(dm_xh, "abx2")
344 call this%aby2%init(dm_xh, "aby2")
345 call this%abz2%init(dm_xh, "abz2")
346 call this%advx%init(dm_xh, "advx")
347 call this%advy%init(dm_xh, "advy")
348 call this%advz%init(dm_xh, "advz")
349 end associate
350
351 call this%du%init(this%dm_Xh, 'du')
352 call this%dv%init(this%dm_Xh, 'dv')
353 call this%dw%init(this%dm_Xh, 'dw')
354 call this%dp%init(this%dm_Xh, 'dp')
355
356 ! Set up boundary conditions
357 call this%setup_bcs(user, params)
358
359 ! Check if we need to output boundaries
360 call json_get_or_default(params, 'case.output_boundary', found, .false.)
361 if (found) call this%write_boundary_conditions()
362
363 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
364 this%pr_projection_activ_step)
365
366 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
367 this%vel_projection_activ_step)
368
369
370
371 ! Determine the time-interpolation scheme
372 call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
373 if (params%valid_path('case.fluid.flow_rate_force')) then
374 call this%vol_flow%init(this%dm_Xh, params)
375 end if
376
377 ! Setup pressure solver
378 call neko_log%section("Pressure solver")
379
380 call json_get_or_default(params, &
381 'case.fluid.pressure_solver.max_iterations', &
382 solver_maxiter, 800)
383 call json_get(params, 'case.fluid.pressure_solver.type', solver_type)
384 call json_get(params, 'case.fluid.pressure_solver.preconditioner.type', &
385 precon_type)
386 call json_get(params, &
387 'case.fluid.pressure_solver.preconditioner', precon_params)
388 call json_get(params, 'case.fluid.pressure_solver.absolute_tolerance', &
389 abs_tol)
390 call json_get_or_default(params, 'case.fluid.pressure_solver.monitor', &
391 monitor, .false.)
392 call neko_log%message('Type : ('// trim(solver_type) // &
393 ', ' // trim(precon_type) // ')')
394 write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
395 call neko_log%message(log_buf)
396
397 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
398 solver_type, solver_maxiter, abs_tol, monitor)
399 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
400 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
401 precon_type, precon_params)
402 call neko_log%end_section()
403
404 ! Initialize the advection factory
405 call json_get_or_default(params, 'case.fluid.advection', advection, .true.)
406 call json_get(params, 'case.numerics', numerics_params)
407 call advection_factory(this%adv, numerics_params, this%c_Xh, &
408 this%ulag, this%vlag, this%wlag, &
409 chkp%dtlag, chkp%tlag, this%ext_bdf, &
410 .not. advection)
411 ! Should be in init_base maybe?
412 this%chkp => chkp
413 ! This is probably scheme specific
414 ! Should not be init really, but more like, add fluid or something...
415 call this%chkp%init(this%u, this%v, this%w, this%p)
416
417 this%chkp%abx1 => this%abx1
418 this%chkp%abx2 => this%abx2
419 this%chkp%aby1 => this%aby1
420 this%chkp%aby2 => this%aby2
421 this%chkp%abz1 => this%abz1
422 this%chkp%abz2 => this%abz2
423 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
424
425
426
427 call neko_log%end_section()
428
429 end subroutine fluid_pnpn_init
430
431 subroutine fluid_pnpn_restart(this, chkp)
432 class(fluid_pnpn_t), target, intent(inout) :: this
433 type(chkp_t), intent(inout) :: chkp
434 real(kind=rp) :: dtlag(10), tlag(10)
435 type(field_t) :: u_temp, v_temp, w_temp
436 integer :: i, j, n
437
438 dtlag = chkp%dtlag
439 tlag = chkp%tlag
440
441 n = this%u%dof%size()
442 if (allocated(chkp%previous_mesh%elements) .or. &
443 chkp%previous_Xh%lx .ne. this%Xh%lx) then
444 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
445 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
446 wlag => this%wlag)
447 do concurrent(j = 1:n)
448 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
449 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
450 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
451 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
452 end do
453 do i = 1, this%ulag%size()
454 do concurrent(j = 1:n)
455 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
456 * c_xh%mult(j,1,1,1)
457 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
458 * c_xh%mult(j,1,1,1)
459 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
460 * c_xh%mult(j,1,1,1)
461 end do
462 end do
463 end associate
464 end if
465
466 if (neko_bcknd_device .eq. 1) then
467 associate(u => this%u, v => this%v, w => this%w, &
468 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
469 p => this%p)
470 call device_memcpy(u%x, u%x_d, u%dof%size(), &
471 host_to_device, sync = .false.)
472 call device_memcpy(v%x, v%x_d, v%dof%size(), &
473 host_to_device, sync = .false.)
474 call device_memcpy(w%x, w%x_d, w%dof%size(), &
475 host_to_device, sync = .false.)
476 call device_memcpy(p%x, p%x_d, p%dof%size(), &
477 host_to_device, sync = .false.)
478 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
479 u%dof%size(), host_to_device, sync = .false.)
480 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
481 u%dof%size(), host_to_device, sync = .false.)
482
483 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
484 v%dof%size(), host_to_device, sync = .false.)
485 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
486 v%dof%size(), host_to_device, sync = .false.)
487
488 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
489 w%dof%size(), host_to_device, sync = .false.)
490 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
491 w%dof%size(), host_to_device, sync = .false.)
492 call device_memcpy(this%abx1%x, this%abx1%x_d, &
493 w%dof%size(), host_to_device, sync = .false.)
494 call device_memcpy(this%abx2%x, this%abx2%x_d, &
495 w%dof%size(), host_to_device, sync = .false.)
496 call device_memcpy(this%aby1%x, this%aby1%x_d, &
497 w%dof%size(), host_to_device, sync = .false.)
498 call device_memcpy(this%aby2%x, this%aby2%x_d, &
499 w%dof%size(), host_to_device, sync = .false.)
500 call device_memcpy(this%abz1%x, this%abz1%x_d, &
501 w%dof%size(), host_to_device, sync = .false.)
502 call device_memcpy(this%abz2%x, this%abz2%x_d, &
503 w%dof%size(), host_to_device, sync = .false.)
504 call device_memcpy(this%advx%x, this%advx%x_d, &
505 w%dof%size(), host_to_device, sync = .false.)
506 call device_memcpy(this%advy%x, this%advy%x_d, &
507 w%dof%size(), host_to_device, sync = .false.)
508 call device_memcpy(this%advz%x, this%advz%x_d, &
509 w%dof%size(), host_to_device, sync = .false.)
510 end associate
511 end if
512 ! Make sure that continuity is maintained (important for interpolation)
513 ! Do not do this for lagged rhs
514 ! (derivatives are not necessairly coninous across elements)
515
516 if (allocated(chkp%previous_mesh%elements) &
517 .or. chkp%previous_Xh%lx .ne. this%Xh%lx) then
518 call this%gs_Xh%op(this%u, gs_op_add)
519 call this%gs_Xh%op(this%v, gs_op_add)
520 call this%gs_Xh%op(this%w, gs_op_add)
521 call this%gs_Xh%op(this%p, gs_op_add)
522
523 do i = 1, this%ulag%size()
524 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
525 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
526 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
527 end do
528 end if
529
530 end subroutine fluid_pnpn_restart
531
532 subroutine fluid_pnpn_free(this)
533 class(fluid_pnpn_t), intent(inout) :: this
534
535 !Deallocate velocity and pressure fields
536 call this%scheme_free()
537
538 if (allocated(this%ext_bdf)) then
539 call this%ext_bdf%free()
540 deallocate(this%ext_bdf)
541 end if
542
543 call this%bc_vel_res%free()
544 call this%bc_du%free()
545 call this%bc_dv%free()
546 call this%bc_dw%free()
547 call this%bc_dp%free()
548
549 call this%bc_prs_surface%free()
550 call this%bc_sym_surface%free()
551 call this%bclst_vel_res%free()
552 call this%bclst_dp%free()
553 call this%proj_prs%free()
554 call this%proj_vel%free()
555
556 call this%p_res%free()
557 call this%u_res%free()
558 call this%v_res%free()
559 call this%w_res%free()
560
561 call this%du%free()
562 call this%dv%free()
563 call this%dw%free()
564 call this%dp%free()
565
566 call this%abx1%free()
567 call this%aby1%free()
568 call this%abz1%free()
569
570 call this%abx2%free()
571 call this%aby2%free()
572 call this%abz2%free()
573
574 call this%advx%free()
575 call this%advy%free()
576 call this%advz%free()
577
578 if (allocated(this%adv)) then
579 call this%adv%free()
580 deallocate(this%adv)
581 end if
582
583 if (allocated(this%Ax_vel)) then
584 deallocate(this%Ax_vel)
585 end if
586
587 if (allocated(this%Ax_prs)) then
588 deallocate(this%Ax_prs)
589 end if
590
591 if (allocated(this%prs_res)) then
592 deallocate(this%prs_res)
593 end if
594
595 if (allocated(this%vel_res)) then
596 deallocate(this%vel_res)
597 end if
598
599 if (allocated(this%sumab)) then
600 deallocate(this%sumab)
601 end if
602
603 if (allocated(this%makeabf)) then
604 deallocate(this%makeabf)
605 end if
606
607 if (allocated(this%makebdf)) then
608 deallocate(this%makebdf)
609 end if
610
611 if (allocated(this%makeoifs)) then
612 deallocate(this%makeoifs)
613 end if
614
615 if (allocated(this%ext_bdf)) then
616 deallocate(this%ext_bdf)
617 end if
618
619 call this%vol_flow%free()
620
621 end subroutine fluid_pnpn_free
622
629 subroutine fluid_pnpn_step(this, time, dt_controller)
630 class(fluid_pnpn_t), target, intent(inout) :: this
631 type(time_state_t), intent(in) :: time
632 type(time_step_controller_t), intent(in) :: dt_controller
633 ! number of degrees of freedom
634 integer :: n
635 ! Solver results monitors (pressure + 3 velocity)
636 type(ksp_monitor_t) :: ksp_results(4)
637
638 type(file_t) :: dump_file
639 class(bc_t), pointer :: bc_i
640 type(non_normal_t), pointer :: bc_j
641
642 if (this%freeze) return
643
644 n = this%dm_Xh%size()
645
646 call profiler_start_region('Fluid', 1)
647 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
648 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
649 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
650 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
651 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
652 xh => this%Xh, &
653 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
654 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
655 msh => this%msh, prs_res => this%prs_res, &
656 source_term => this%source_term, vel_res => this%vel_res, &
657 sumab => this%sumab, makeoifs => this%makeoifs, &
658 makeabf => this%makeabf, makebdf => this%makebdf, &
659 vel_projection_dim => this%vel_projection_dim, &
660 pr_projection_dim => this%pr_projection_dim, &
661 oifs => this%oifs, &
662 rho => this%rho, mu_tot => this%mu_tot, &
663 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
664 t => time%t, tstep => time%tstep, dt => time%dt, &
665 ext_bdf => this%ext_bdf, event => glb_cmd_event)
666
667 ! Extrapolate the velocity if it's not done in nut_field estimation
668 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
669 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
670
671 ! Compute the source terms
672 call this%source_term%compute(time)
673
674 ! Add Neumann bc contributions to the RHS
675 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
676 this%dm_Xh%size(), time, strong = .false.)
677
678 if (oifs) then
679 ! Add the advection operators to the right-hand-side.
680 call this%adv%compute(u, v, w, &
681 this%advx, this%advy, this%advz, &
682 xh, this%c_Xh, dm_xh%size(), dt)
683
684 ! At this point the RHS contains the sum of the advection operator and
685 ! additional source terms, evaluated using the velocity field from the
686 ! previous time-step. Now, this value is used in the explicit time
687 ! scheme to advance both terms in time.
688
689 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
690 this%abx2, this%aby2, this%abz2, &
691 f_x%x, f_y%x, f_z%x, &
692 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
693
694 ! Now, the source terms from the previous time step are added to the RHS.
695 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
696 f_x%x, f_y%x, f_z%x, &
697 rho%x(1,1,1,1), dt, n)
698 else
699 ! Add the advection operators to the right-hand-side.
700 call this%adv%compute(u, v, w, &
701 f_x, f_y, f_z, &
702 xh, this%c_Xh, dm_xh%size())
703
704 ! At this point the RHS contains the sum of the advection operator and
705 ! additional source terms, evaluated using the velocity field from the
706 ! previous time-step. Now, this value is used in the explicit time
707 ! scheme to advance both terms in time.
708 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
709 this%abx2, this%aby2, this%abz2, &
710 f_x%x, f_y%x, f_z%x, &
711 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
712
713 ! Add the RHS contributions coming from the BDF scheme.
714 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
715 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
716 ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
717 end if
718
719 call ulag%update()
720 call vlag%update()
721 call wlag%update()
722
723 call this%bc_apply_vel(time, strong = .true.)
724 call this%bc_apply_prs(time)
725
726 ! Update material properties if necessary
727 call this%update_material_properties(time)
728
729 ! Compute pressure residual.
730 call profiler_start_region('Pressure_residual', 18)
731
732 call prs_res%compute(p, p_res,&
733 u, v, w, &
734 u_e, v_e, w_e, &
735 f_x, f_y, f_z, &
736 c_xh, gs_xh, &
737 this%bc_prs_surface, this%bc_sym_surface,&
738 ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
739 mu_tot, rho, event)
740
741 ! De-mean the pressure residual when no strong pressure boundaries present
742 if (.not. this%prs_dirichlet) call ortho(p_res%x, this%glb_n_points, n)
743
744 call gs_xh%op(p_res, gs_op_add, event)
745 call device_event_sync(event)
746
747 ! Set the residual to zero at strong pressure boundaries.
748 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
749
750
751 call profiler_end_region('Pressure_residual', 18)
752
753
754 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
755 'Pressure')
756
757 call this%pc_prs%update()
758
759 call profiler_start_region('Pressure_solve', 3)
760
761 ! Solve for the pressure increment.
762 ksp_results(1) = &
763 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
764 this%bclst_dp, gs_xh)
765 ksp_results(1)%name = 'Pressure'
766
767
768 call profiler_end_region('Pressure_solve', 3)
769
770 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
771 this%bclst_dp, gs_xh, n, tstep, dt_controller)
772
773 ! Update the pressure with the increment. Demean if necessary.
774 call field_add2(p, dp, n)
775 if (.not. this%prs_dirichlet) call ortho(p%x, this%glb_n_points, n)
776
777 ! Compute velocity residual.
778 call profiler_start_region('Velocity_residual', 19)
779 call vel_res%compute(ax_vel, u, v, w, &
780 u_res, v_res, w_res, &
781 p, &
782 f_x, f_y, f_z, &
783 c_xh, msh, xh, &
784 mu_tot, rho, ext_bdf%diffusion_coeffs(1), &
785 dt, dm_xh%size())
786
787 call gs_xh%op(u_res, gs_op_add, event)
788 call device_event_sync(event)
789 call gs_xh%op(v_res, gs_op_add, event)
790 call device_event_sync(event)
791 call gs_xh%op(w_res, gs_op_add, event)
792 call device_event_sync(event)
793
794 ! Set residual to zero at strong velocity boundaries.
795 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
796
797
798 call profiler_end_region('Velocity_residual', 19)
799
800 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
801 tstep, c_xh, n, dt_controller, 'Velocity')
802
803 call this%pc_vel%update()
804
805 call profiler_start_region("Velocity_solve", 4)
806 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
807 u_res%x, v_res%x, w_res%x, n, c_xh, &
808 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
809 this%ksp_vel%max_iter)
810 call profiler_end_region("Velocity_solve", 4)
811 if (this%full_stress_formulation) then
812 ksp_results(2)%name = 'Momentum'
813 else
814 ksp_results(2)%name = 'X-Velocity'
815 ksp_results(3)%name = 'Y-Velocity'
816 ksp_results(4)%name = 'Z-Velocity'
817 end if
818
819 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
820 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
821 dt_controller)
822
823 if (neko_bcknd_device .eq. 1) then
824 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
825 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
826 else
827 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
828 end if
829
830 if (this%forced_flow_rate) then
831 ! Horrible mu hack?!
832 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
833 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
834 dt, this%bclst_dp, this%bclst_du, this%bclst_dv, &
835 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
836 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
837 this%ksp_vel%max_iter)
838 end if
839
840 call fluid_step_info(time, ksp_results, &
841 this%full_stress_formulation, this%strict_convergence)
842
843 end associate
844 call profiler_end_region('Fluid', 1)
845 end subroutine fluid_pnpn_step
846
849 subroutine fluid_pnpn_setup_bcs(this, user, params)
850 class(fluid_pnpn_t), target, intent(inout) :: this
851 type(user_t), target, intent(in) :: user
852 type(json_file), intent(inout) :: params
853 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
854 class(bc_t), pointer :: bc_i
855 type(json_core) :: core
856 type(json_value), pointer :: bc_object
857 type(json_file) :: bc_subdict
858 logical :: found
859 ! Monitor which boundary zones have been marked
860 logical, allocatable :: marked_zones(:)
861 integer, allocatable :: zone_indices(:)
862
863 ! Lists for the residuals and solution increments
864 call this%bclst_vel_res%init()
865 call this%bclst_du%init()
866 call this%bclst_dv%init()
867 call this%bclst_dw%init()
868 call this%bclst_dp%init()
869
870 call this%bc_vel_res%init_from_components(this%c_Xh)
871 call this%bc_du%init_from_components(this%c_Xh)
872 call this%bc_dv%init_from_components(this%c_Xh)
873 call this%bc_dw%init_from_components(this%c_Xh)
874 call this%bc_dp%init_from_components(this%c_Xh)
875
876 ! Special PnPn boundary conditions for pressure
877 call this%bc_prs_surface%init_from_components(this%c_Xh)
878 call this%bc_sym_surface%init_from_components(this%c_Xh)
879
880 ! Populate bcs_vel and bcs_prs based on the case file
881 if (params%valid_path('case.fluid.boundary_conditions')) then
882 call params%info('case.fluid.boundary_conditions', n_children = n_bcs)
883 call params%get_core(core)
884 call params%get('case.fluid.boundary_conditions', bc_object, found)
885
886 !
887 ! Velocity bcs
888 !
889 call this%bcs_vel%init(n_bcs)
890
891 allocate(marked_zones(size(this%msh%labeled_zones)))
892 marked_zones = .false.
893
894 do i = 1, n_bcs
895 ! Create a new json containing just the subdict for this bc
896 call json_extract_item(core, bc_object, i, bc_subdict)
897
898 call json_get(bc_subdict, "zone_indices", zone_indices)
899
900 ! Check that we are not trying to assing a bc to zone, for which one
901 ! has already been assigned and that the zone has more than 0 size
902 ! in the mesh.
903 do j = 1, size(zone_indices)
904 zone_size = this%msh%labeled_zones(zone_indices(j))%size
905 call mpi_allreduce(zone_size, global_zone_size, 1, &
906 mpi_integer, mpi_max, neko_comm, ierr)
907
908 if (global_zone_size .eq. 0) then
909 write(error_unit, '(A, A, I0, A, A, I0, A)') "*** ERROR ***: ",&
910 "Zone index ", zone_indices(j), &
911 " is invalid as this zone has 0 size, meaning it ", &
912 "is not in the mesh. Check fluid boundary condition ", &
913 i, "."
914 error stop
915 end if
916
917 if (marked_zones(zone_indices(j)) .eqv. .true.) then
918 write(error_unit, '(A, A, I0, A, A, A, A)') "*** ERROR ***: ", &
919 "Zone with index ", zone_indices(j), &
920 " has already been assigned a boundary condition. ", &
921 "Please check your boundary_conditions entry for the ", &
922 "fluid and make sure that each zone index appears only ", &
923 "in a single boundary condition."
924 error stop
925 else
926 marked_zones(zone_indices(j)) = .true.
927 end if
928 end do
929
930 bc_i => null()
931 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
932
933 ! Not all bcs require an allocation for velocity in particular,
934 ! so we check.
935 if (associated(bc_i)) then
936
937 ! We need to treat mixed bcs separately because they are by
938 ! convention marked weak and currently contain nested
939 ! bcs, some of which are strong.
940 select type (bc_i)
941 type is (symmetry_t)
942 ! Symmetry has 3 internal bcs, but only one actually contains
943 ! markings.
944 ! Symmetry's apply_scalar doesn't do anything, so we need to mark
945 ! individual nested bcs to the du,dv,dw, whereas the vel_res can
946 ! just get symmetry as a whole, because on this list we call
947 ! apply_vector.
948 ! Additionally we have to mark the special surface bc for p.
949 call this%bclst_vel_res%append(bc_i)
950 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
951 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
952 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
953
954 call this%bcs_vel%append(bc_i)
955
956 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
957 type is (non_normal_t)
958 ! This is a bc for the residuals and increments, not the
959 ! velocity itself. So, don't append to bcs_vel
960 call this%bclst_vel_res%append(bc_i)
961 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
962 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
963 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
964 type is (shear_stress_t)
965 ! Same as symmetry
966 call this%bclst_vel_res%append(bc_i%symmetry)
967 call this%bclst_du%append(bc_i%symmetry%bc_x)
968 call this%bclst_dv%append(bc_i%symmetry%bc_y)
969 call this%bclst_dw%append(bc_i%symmetry%bc_z)
970
971 call this%bcs_vel%append(bc_i)
972 type is (wall_model_bc_t)
973 ! Same as symmetry
974 call this%bclst_vel_res%append(bc_i%symmetry)
975 call this%bclst_du%append(bc_i%symmetry%bc_x)
976 call this%bclst_dv%append(bc_i%symmetry%bc_y)
977 call this%bclst_dw%append(bc_i%symmetry%bc_z)
978
979 call this%bcs_vel%append(bc_i)
980 class default
981
982 ! For the default case we use our dummy zero_dirichlet bcs to
983 ! mark the same faces as in ordinary velocity dirichlet
984 ! conditions.
985 ! Additionally we mark the special PnPn pressure bc.
986 if (bc_i%strong .eqv. .true.) then
987 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
988 call this%bc_du%mark_facets(bc_i%marked_facet)
989 call this%bc_dv%mark_facets(bc_i%marked_facet)
990 call this%bc_dw%mark_facets(bc_i%marked_facet)
991
992 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
993 end if
994
995 call this%bcs_vel%append(bc_i)
996 end select
997 end if
998 end do
999
1000 ! Make sure all labeled zones with non-zero size have been marked
1001 do i = 1, size(this%msh%labeled_zones)
1002 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1003 (marked_zones(i) .eqv. .false.)) then
1004 write(error_unit, '(A, A, I0)') "*** ERROR ***: ", &
1005 "No fluid boundary condition assigned to zone ", i
1006 error stop
1007 end if
1008 end do
1009
1010 !
1011 ! Pressure bcs
1012 !
1013 call this%bcs_prs%init(n_bcs)
1014
1015 do i = 1, n_bcs
1016 ! Create a new json containing just the subdict for this bc
1017 call json_extract_item(core, bc_object, i, bc_subdict)
1018 bc_i => null()
1019 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1020
1021 ! Not all bcs require an allocation for pressure in particular,
1022 ! so we check.
1023 if (associated(bc_i)) then
1024 call this%bcs_prs%append(bc_i)
1025
1026 ! Mark strong bcs in the dummy dp bc to force zero change.
1027 if (bc_i%strong .eqv. .true.) then
1028 call this%bc_dp%mark_facets(bc_i%marked_facet)
1029 end if
1030
1031 end if
1032
1033 end do
1034 else
1035 ! Check that there are no labeled zones, i.e. all are periodic.
1036 do i = 1, size(this%msh%labeled_zones)
1037 if (this%msh%labeled_zones(i)%size .gt. 0) then
1038 call neko_error("No boundary_conditions entry in the case file!")
1039 end if
1040 end do
1041
1042 ! For a pure periodic case, we still need to initilise the bc lists
1043 ! to a zero size to avoid issues with apply() in step()
1044 call this%bcs_vel%init()
1045 call this%bcs_prs%init()
1046
1047 end if
1048
1049 call this%bc_prs_surface%finalize()
1050 call this%bc_sym_surface%finalize()
1051
1052 call this%bc_vel_res%finalize()
1053 call this%bc_du%finalize()
1054 call this%bc_dv%finalize()
1055 call this%bc_dw%finalize()
1056 call this%bc_dp%finalize()
1057
1058 call this%bclst_vel_res%append(this%bc_vel_res)
1059 call this%bclst_du%append(this%bc_du)
1060 call this%bclst_dv%append(this%bc_dv)
1061 call this%bclst_dw%append(this%bc_dw)
1062 call this%bclst_dp%append(this%bc_dp)
1063
1064 ! If we have no strong pressure bcs, we will demean the pressure
1065 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1066 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1067 mpi_logical, mpi_lor, neko_comm)
1068
1069 end subroutine fluid_pnpn_setup_bcs
1070
1072 subroutine fluid_pnpn_write_boundary_conditions(this)
1073 use inflow, only : inflow_t
1075 use blasius, only : blasius_t
1077 use dong_outflow, only : dong_outflow_t
1078 class(fluid_pnpn_t), target, intent(inout) :: this
1079 type(dirichlet_t) :: bdry_mask
1080 type(field_t), pointer :: bdry_field
1081 type(file_t) :: bdry_file
1082 integer :: temp_index, i
1083 class(bc_t), pointer :: bci
1084 character(len=LOG_SIZE) :: log_buf
1085
1086 call neko_log%section("Fluid boundary conditions")
1087 write(log_buf, '(A)') 'Marking using integer keys in bdry0.f00000'
1088 call neko_log%message(log_buf)
1089 write(log_buf, '(A)') 'Condition-value pairs: '
1090 call neko_log%message(log_buf)
1091 write(log_buf, '(A)') ' no_slip = 1'
1092 call neko_log%message(log_buf)
1093 write(log_buf, '(A)') ' velocity_value = 2'
1094 call neko_log%message(log_buf)
1095 write(log_buf, '(A)') ' outflow, normal_outflow (+dong) = 3'
1096 call neko_log%message(log_buf)
1097 write(log_buf, '(A)') ' symmetry = 4'
1098 call neko_log%message(log_buf)
1099 write(log_buf, '(A)') ' periodic = 6'
1100 call neko_log%message(log_buf)
1101 write(log_buf, '(A)') ' user_velocity = 7'
1102 call neko_log%message(log_buf)
1103 write(log_buf, '(A)') ' user_pressure = 8'
1104 call neko_log%message(log_buf)
1105 write(log_buf, '(A)') ' shear_stress = 9'
1106 call neko_log%message(log_buf)
1107 write(log_buf, '(A)') ' wall_modelling = 10'
1108 call neko_log%message(log_buf)
1109 write(log_buf, '(A)') ' blasius_profile = 11'
1110 call neko_log%message(log_buf)
1111 call neko_log%end_section()
1112
1113 call this%scratch%request_field(bdry_field, temp_index)
1114 bdry_field = 0.0_rp
1115
1116
1117
1118 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1119 call bdry_mask%mark_zone(this%msh%periodic)
1120 call bdry_mask%finalize()
1121 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1122 call bdry_mask%free()
1123
1124 do i = 1, this%bcs_prs%size()
1125 bci => this%bcs_prs%get(i)
1126 select type (bc => bci)
1127 type is (zero_dirichlet_t)
1128 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1129 call bdry_mask%mark_facets(bci%marked_facet)
1130 call bdry_mask%finalize()
1131 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1132 call bdry_mask%free()
1133 type is (dong_outflow_t)
1134 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1135 call bdry_mask%mark_facets(bci%marked_facet)
1136 call bdry_mask%finalize()
1137 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1138 call bdry_mask%free()
1139 type is (field_dirichlet_t)
1140 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1141 call bdry_mask%mark_facets(bci%marked_facet)
1142 call bdry_mask%finalize()
1143 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1144 call bdry_mask%free()
1145 end select
1146 end do
1147
1148 do i = 1, this%bcs_vel%size()
1149 bci => this%bcs_vel%get(i)
1150 select type (bc => bci)
1151 type is (zero_dirichlet_t)
1152 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1153 call bdry_mask%mark_facets(bci%marked_facet)
1154 call bdry_mask%finalize()
1155 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1156 call bdry_mask%free()
1157 type is (inflow_t)
1158 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1159 call bdry_mask%mark_facets(bci%marked_facet)
1160 call bdry_mask%finalize()
1161 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1162 call bdry_mask%free()
1163 type is (symmetry_t)
1164 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1165 call bdry_mask%mark_facets(bci%marked_facet)
1166 call bdry_mask%finalize()
1167 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1168 call bdry_mask%free()
1169 type is (field_dirichlet_vector_t)
1170 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1171 call bdry_mask%mark_facets(bci%marked_facet)
1172 call bdry_mask%finalize()
1173 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1174 call bdry_mask%free()
1175 type is (shear_stress_t)
1176 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1177 call bdry_mask%mark_facets(bci%marked_facet)
1178 call bdry_mask%finalize()
1179 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1180 call bdry_mask%free()
1181 type is (wall_model_bc_t)
1182 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1183 call bdry_mask%mark_facets(bci%marked_facet)
1184 call bdry_mask%finalize()
1185 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1186 call bdry_mask%free()
1187 type is (blasius_t)
1188 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1189 call bdry_mask%mark_facets(bci%marked_facet)
1190 call bdry_mask%finalize()
1191 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1192 call bdry_mask%free()
1193 end select
1194 end do
1195
1196
1197 call bdry_file%init('bdry.fld')
1198 call bdry_file%write(bdry_field)
1199
1200 call this%scratch%relinquish_field(temp_index)
1201 end subroutine fluid_pnpn_write_boundary_conditions
1202
1203end module fluid_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 Blasius profile dirichlet condition.
Definition blasius.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_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
subroutine, public device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
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
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a dong outflow condition.
Dirichlet condition applied in the facet normal direction.
Defines inflow dirichlet conditions.
Defines user dirichlet condition for a scalar field.
subroutine, public field_add2(a, b, n)
Vector addition .
subroutine, public field_copy(a, b, n)
Copy a vector .
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Auxiliary routines for fluid solvers.
Definition fluid_aux.f90:2
subroutine, public fluid_step_info(time, ksp_results, full_stress_formulation, strict_convergence)
Prints for prs, velx, vely, velz the following: Number of iterations, start residual,...
Definition fluid_aux.f90:19
Modular version of the Classic Nek5000 Pn/Pn formulation for fluids.
subroutine fluid_pnpn_setup_bcs(this, user, params)
Sets up the boundary condition for the scheme.
subroutine fluid_pnpn_restart(this, chkp)
subroutine fluid_pnpn_write_boundary_conditions(this)
Write a field with boundary condition specifications.
subroutine fluid_pnpn_step(this, time, dt_controller)
Advance fluid simulation in time.
subroutine fluid_pnpn_free(this)
subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
Boundary condition factory for pressure.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Defines inflow dirichlet conditions.
Definition inflow.f90:34
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:70
integer, parameter, public log_size
Definition log.f90:46
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition mathops.f90:65
subroutine, public opcolv(a1, a2, a3, c, gdim, n)
Definition mathops.f90:97
subroutine, public opadd2cm(a1, a2, a3, b1, b2, b3, c, n, gdim)
Definition mathops.f90:142
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
Dirichlet condition on axis aligned plane in the non normal direction.
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public ortho(x, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
Defines Pressure and velocity residuals in the Pn-Pn formulation.
Definition pnpn_res.f90:34
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 Couple projections for velocity.
Project x onto X, the space of old solutions and back again.
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition rhs_maker.f90:38
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_type_error(base_type, wrong_type, known_types)
Reports an error allocating a type for a particular base pointer class.
Definition utils.f90:313
Defines the wall_model_bc_t type. Maintainer: Timofey Mukha.
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
Blasius profile for inlet (vector valued).
Definition blasius.f90:53
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:48
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...
Dirichlet condition in facet normal direction.
User defined dirichlet condition, for which the user can work with an entire field....
Extension of the user defined dirichlet condition field_dirichlet
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:55
Dirichlet condition for inlet (vector valued)
Definition inflow.f90:47
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Dirichlet condition in non normal direction of a plane.
Abstract type to compute pressure residual.
Definition pnpn_res.f90:48
Abstract type to compute velocity residual.
Definition pnpn_res.f90:54
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 extrapolated velocity field for the pressure equation.
Definition rhs_maker.f90:46
A shear stress boundary condition.
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:50
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...
A shear stress boundary condition, computing the stress values using a wall model.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...