Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
zero_dirichlet.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
37 use num_types, only : rp
38 use bc, only : bc_t
39 use, intrinsic :: iso_c_binding, only : c_ptr
40 use coefs, only : coef_t
41 use json_module, only : json_file
42 use time_state, only : time_state_t
43 implicit none
44 private
45
49 type, public, extends(bc_t) :: zero_dirichlet_t
50 contains
51 procedure, pass(this) :: apply_scalar => zero_dirichlet_apply_scalar
52 procedure, pass(this) :: apply_vector => zero_dirichlet_apply_vector
53 procedure, pass(this) :: apply_scalar_dev => &
55 procedure, pass(this) :: apply_vector_dev => &
58 procedure, pass(this) :: init => zero_dirichlet_init
60 procedure, pass(this) :: init_from_components => &
63 procedure, pass(this) :: free => zero_dirichlet_free
65 procedure, pass(this) :: finalize => zero_dirichlet_finalize
66 end type zero_dirichlet_t
67
68contains
69
73 subroutine zero_dirichlet_init(this, coef, json)
74 class(zero_dirichlet_t), intent(inout), target :: this
75 type(coef_t), target, intent(in) :: coef
76 type(json_file), intent(inout) :: json
77
78 call this%init_from_components(coef)
79 end subroutine zero_dirichlet_init
80
84 class(zero_dirichlet_t), intent(inout), target :: this
85 type(coef_t), target, intent(in) :: coef
86
87 call this%init_base(coef)
89
92 subroutine zero_dirichlet_apply_scalar(this, x, n, time, strong)
93 class(zero_dirichlet_t), intent(inout) :: this
94 integer, intent(in) :: n
95 real(kind=rp), intent(inout), dimension(n) :: x
96 type(time_state_t), intent(in), optional :: time
97 logical, intent(in), optional :: strong
98 integer :: i, m, k
99 logical :: strong_
100
101 if (present(strong)) then
102 strong_ = strong
103 else
104 strong_ = .true.
105 end if
106
107 m = this%msk(0)
108
109 if (strong_) then
110 do i = 1, m
111 k = this%msk(i)
112 x(k) = 0d0
113 end do
114 end if
115
116 end subroutine zero_dirichlet_apply_scalar
117
119 subroutine zero_dirichlet_apply_vector(this, x, y, z, n, time, strong)
120 class(zero_dirichlet_t), intent(inout) :: this
121 integer, intent(in) :: n
122 real(kind=rp), intent(inout), dimension(n) :: x
123 real(kind=rp), intent(inout), dimension(n) :: y
124 real(kind=rp), intent(inout), dimension(n) :: z
125 type(time_state_t), intent(in), optional :: time
126 logical, intent(in), optional :: strong
127 integer :: i, m, k
128 logical :: strong_
129
130 if (present(strong)) then
131 strong_ = strong
132 else
133 strong_ = .true.
134 end if
135
136 if (strong_) then
137 m = this%msk(0)
138 do i = 1, m
139 k = this%msk(i)
140 x(k) = 0d0
141 y(k) = 0d0
142 z(k) = 0d0
143 end do
144 end if
145
146 end subroutine zero_dirichlet_apply_vector
147
149 subroutine zero_dirichlet_apply_scalar_dev(this, x_d, time, strong, strm)
150 class(zero_dirichlet_t), intent(inout), target :: this
151 type(c_ptr), intent(inout) :: x_d
152 type(time_state_t), intent(in), optional :: time
153 logical, intent(in), optional :: strong
154 type(c_ptr), intent(inout) :: strm
155 logical :: strong_
156
157 if (present(strong)) then
158 strong_ = strong
159 else
160 strong_ = .true.
161 end if
162
163 if (strong_ .and. (this%msk(0) .gt. 0)) then
164 call device_zero_dirichlet_apply_scalar(this%msk_d, x_d, &
165 size(this%msk), strm)
166 end if
167
169
171 subroutine zero_dirichlet_apply_vector_dev(this, x_d, y_d, z_d, time, &
172 strong, strm)
173 class(zero_dirichlet_t), intent(inout), target :: this
174 type(c_ptr), intent(inout) :: x_d
175 type(c_ptr), intent(inout) :: y_d
176 type(c_ptr), intent(inout) :: z_d
177 type(time_state_t), intent(in), optional :: time
178 logical, intent(in), optional :: strong
179 type(c_ptr), intent(inout) :: strm
180 logical :: strong_
181
182 if (present(strong)) then
183 strong_ = strong
184 else
185 strong_ = .true.
186 end if
187
188 if (strong_ .and. (this%msk(0) .gt. 0)) then
189 call device_zero_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, &
190 size(this%msk), strm)
191 end if
192
194
196 subroutine zero_dirichlet_free(this)
197 class(zero_dirichlet_t), target, intent(inout) :: this
198
199 call this%free_base()
200
201 end subroutine zero_dirichlet_free
202
204 subroutine zero_dirichlet_finalize(this, only_facets)
205 class(zero_dirichlet_t), target, intent(inout) :: this
206 logical, optional, intent(in) :: only_facets
207 logical :: only_facets_
208
209 if (present(only_facets)) then
210 only_facets_ = only_facets
211 else
212 only_facets_ = .false.
213 end if
214
215 call this%finalize_base(only_facets_)
216 end subroutine zero_dirichlet_finalize
217
218end module zero_dirichlet
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_zero_dirichlet_apply_vector(msk, x, y, z, m, strm)
subroutine, public device_zero_dirichlet_apply_scalar(msk, x, m, strm)
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Defines a zero-valued Dirichlet boundary condition.
subroutine zero_dirichlet_apply_scalar(this, x, n, time, strong)
Apply boundary condition to a scalar field. to a vector x.
subroutine zero_dirichlet_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Apply boundary condition to a vector field, device version.
subroutine zero_dirichlet_init(this, coef, json)
Constructor.
subroutine zero_dirichlet_apply_scalar_dev(this, x_d, time, strong, strm)
Apply boundary condition to a scalar field, device version.
subroutine zero_dirichlet_free(this)
Destructor.
subroutine zero_dirichlet_init_from_components(this, coef)
Constructor.
subroutine zero_dirichlet_finalize(this, only_facets)
Finalize.
subroutine zero_dirichlet_apply_vector(this, x, y, z, n, time, strong)
Apply boundary condition to a vector field.
Base type for a boundary condition.
Definition bc.f90:62
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
A struct that contains all info about the time, expand as needed.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...