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