Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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 coefs, only : coef_t
40 use json_module, only : json_file
42 use, intrinsic :: iso_c_binding, only : c_ptr
43 use time_state, only : time_state_t
44 implicit none
45 private
46
49 type, public, extends(bc_t) :: dirichlet_t
50 real(kind=rp), private :: g
51 contains
52 procedure, pass(this) :: apply_scalar => dirichlet_apply_scalar
53 procedure, pass(this) :: apply_vector => dirichlet_apply_vector
54 procedure, pass(this) :: apply_scalar_dev => dirichlet_apply_scalar_dev
55 procedure, pass(this) :: apply_vector_dev => dirichlet_apply_vector_dev
56 procedure, pass(this) :: set_g => dirichlet_set_g
58 procedure, pass(this) :: init => dirichlet_init
60 procedure, pass(this) :: init_from_components => &
63 procedure, pass(this) :: free => dirichlet_free
65 procedure, pass(this) :: finalize => dirichlet_finalize
66 end type dirichlet_t
67
68contains
69
73 subroutine dirichlet_init(this, coef, json)
74 class(dirichlet_t), intent(inout), target :: this
75 type(coef_t), target, intent(in) :: coef
76 type(json_file), intent(inout) ::json
77 real(kind=rp) :: g
78
79 call this%init_base(coef)
80 call json_get_or_lookup(json , "value", g)
81
82 this%g = g
83 end subroutine dirichlet_init
84
88 subroutine dirichlet_init_from_components(this, coef, g)
89 class(dirichlet_t), intent(inout), target :: this
90 type(coef_t), target, intent(in) :: coef
91 real(kind=rp), intent(in) :: g
92
93 call this%init_base(coef)
94 this%g = g
96
99 subroutine dirichlet_apply_scalar(this, x, n, time, strong)
100 class(dirichlet_t), intent(inout) :: this
101 integer, intent(in) :: n
102 real(kind=rp), intent(inout), dimension(n) :: x
103 type(time_state_t), intent(in), optional :: time
104 logical, intent(in), optional :: strong
105 integer :: i, m, k
106 logical :: strong_
107
108 if (present(strong)) then
109 strong_ = strong
110 else
111 strong_ = .true.
112 end if
113
114 if (strong_) then
115 m = this%msk(0)
116 do i = 1, m
117 k = this%msk(i)
118 x(k) = this%g
119 end do
120 end if
121 end subroutine dirichlet_apply_scalar
122
125 subroutine dirichlet_apply_vector(this, x, y, z, n, time, strong)
126 class(dirichlet_t), intent(inout) :: this
127 integer, intent(in) :: n
128 real(kind=rp), intent(inout), dimension(n) :: x
129 real(kind=rp), intent(inout), dimension(n) :: y
130 real(kind=rp), intent(inout), dimension(n) :: z
131 type(time_state_t), intent(in), optional :: time
132 logical, intent(in), optional :: strong
133 integer :: i, m, k
134 logical :: strong_
135
136 if (present(strong)) then
137 strong_ = strong
138 else
139 strong_ = .true.
140 end if
141
142 if (strong_) then
143 m = this%msk(0)
144 do i = 1, m
145 k = this%msk(i)
146 x(k) = this%g
147 y(k) = this%g
148 z(k) = this%g
149 end do
150 end if
151
152 end subroutine dirichlet_apply_vector
153
156 subroutine dirichlet_apply_scalar_dev(this, x_d, time, strong, strm)
157 class(dirichlet_t), intent(inout), target :: this
158 type(c_ptr), intent(inout) :: x_d
159 type(time_state_t), intent(in), optional :: time
160 logical, intent(in), optional :: strong
161 type(c_ptr), intent(inout) :: strm
162 logical :: strong_
163
164 if (present(strong)) then
165 strong_ = strong
166 else
167 strong_ = .true.
168 end if
169
170 if (strong_ .and. this%msk(0) .gt. 0) then
171 call device_dirichlet_apply_scalar(this%msk_d, x_d, &
172 this%g, size(this%msk), strm)
173 end if
174
175 end subroutine dirichlet_apply_scalar_dev
176
179 subroutine dirichlet_apply_vector_dev(this, x_d, y_d, z_d, &
180 time, strong, strm)
181 class(dirichlet_t), intent(inout), target :: this
182 type(c_ptr), intent(inout) :: x_d
183 type(c_ptr), intent(inout) :: y_d
184 type(c_ptr), intent(inout) :: z_d
185 type(time_state_t), intent(in), optional :: time
186 logical, intent(in), optional :: strong
187 type(c_ptr), intent(inout) :: strm
188 logical :: strong_
189
190 if (present(strong)) then
191 strong_ = strong
192 else
193 strong_ = .true.
194 end if
195
196 if (strong_ .and. this%msk(0) .gt. 0) then
197 call device_dirichlet_apply_vector(this%msk_d, x_d, y_d, z_d, this%g, &
198 size(this%msk), strm)
199 end if
200
201 end subroutine dirichlet_apply_vector_dev
202
204 subroutine dirichlet_set_g(this, g)
205 class(dirichlet_t), intent(inout) :: this
206 real(kind=rp), intent(in) :: g
207
208 this%g = g
209
210 end subroutine dirichlet_set_g
211
213 subroutine dirichlet_free(this)
214 class(dirichlet_t), target, intent(inout) :: this
215
216 call this%free_base
217
218 end subroutine dirichlet_free
219
221 subroutine dirichlet_finalize(this, only_facets)
222 class(dirichlet_t), target, intent(inout) :: this
223 logical, optional, intent(in) :: only_facets
224 logical :: only_facets_
225
226 if (present(only_facets)) then
227 only_facets_ = only_facets
228 else
229 only_facets_ = .false.
230 end if
231
232 call this%finalize_base(only_facets_)
233 end subroutine dirichlet_finalize
234
235end module dirichlet
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_dirichlet_apply_scalar(msk, x, g, m, strm)
subroutine, public device_dirichlet_apply_vector(msk, x, y, z, g, m, strm)
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
subroutine dirichlet_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z (device version)
subroutine dirichlet_finalize(this, only_facets)
Finalize.
subroutine dirichlet_apply_scalar_dev(this, x_d, time, strong, strm)
Boundary condition apply for a generic Dirichlet condition to a vector x (device version)
subroutine dirichlet_apply_scalar(this, x, n, time, strong)
Boundary condition apply for a generic Dirichlet condition to a vector x.
subroutine dirichlet_apply_vector(this, x, y, z, n, time, strong)
Boundary condition apply for a generic Dirichlet condition to vectors x, y and z.
subroutine dirichlet_free(this)
Destructor.
subroutine dirichlet_init_from_components(this, coef, g)
Constructor from components.
Definition dirichlet.f90:89
subroutine dirichlet_init(this, coef, json)
Constructor from JSON.
Definition dirichlet.f90:74
subroutine dirichlet_set_g(this, g)
Set value of .
Utilities for retrieving parameters from the case files.
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
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
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
A struct that contains all info about the time, expand as needed.