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