Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
RAMP_mapping.f90
Go to the documentation of this file.
1! Copyright (c) 2023, 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!
33!
36 use num_types, only: rp
37 use mapping, only: mapping_t
38 use json_module, only: json_file
39 use field, only: field_t
40 use coefs, only: coef_t
41 use neko_config, only: neko_bcknd_device
46 implicit none
47 private
48
76
77 type, public, extends(mapping_t) :: ramp_mapping_t
79 real(kind=rp) :: f_min
81 real(kind=rp) :: f_max
83 real(kind=rp) :: q
86 logical :: convex_up
87
88 contains
90 procedure, pass(this) :: init => ramp_mapping_init_from_json
92 procedure, pass(this) :: init_from_attributes => &
95 procedure, pass(this) :: free => ramp_mapping_free
97 procedure, pass(this) :: apply_forward => ramp_mapping_apply
99 procedure, pass(this) :: apply_backward => &
101 end type ramp_mapping_t
102
103contains
104
106 subroutine ramp_mapping_init_from_json(this, json, coef)
107 class(ramp_mapping_t), intent(inout) :: this
108 type(json_file), intent(inout) :: json
109 type(coef_t), intent(inout) :: coef
110
111 ! do the JSON stuff later
112 this%f_min = 0.0_rp
113 this%f_max = 1000.0_rp
114 this%q = 1.0_rp
115 this%convex_up = .true.
116
117 call this%init_base(json, coef)
118 call ramp_mapping_init_from_attributes(this, coef)
119
120 end subroutine ramp_mapping_init_from_json
121
124 class(ramp_mapping_t), intent(inout) :: this
125 type(coef_t), intent(inout) :: coef
126
127 ! there's actually nothing to do here.
128
130
132 subroutine ramp_mapping_free(this)
133 class(ramp_mapping_t), intent(inout) :: this
134
135 call this%free_base()
136
137 end subroutine ramp_mapping_free
138
142 subroutine ramp_mapping_apply(this, X_out, X_in)
143 class(ramp_mapping_t), intent(inout) :: this
144 type(field_t), intent(in) :: X_in
145 type(field_t), intent(inout) :: X_out
146
147 if (this%convex_up .eqv. .true.) then
148 call convex_up_ramp_mapping_apply(this%f_min, this%f_max, &
149 this%q, x_out, x_in)
150 else
151 call convex_down_ramp_mapping_apply(this%f_min, this%f_max, &
152 this%q, x_out, x_in)
153 end if
154
155 end subroutine ramp_mapping_apply
156
157
162 subroutine ramp_mapping_apply_backward(this, dF_dX_in, dF_dX_out, X_in)
163 class(ramp_mapping_t), intent(inout) :: this
164 type(field_t), intent(in) :: X_in
165 type(field_t), intent(in) :: dF_dX_out
166 type(field_t), intent(inout) :: dF_dX_in
167
168 if (this%convex_up .eqv. .true.) then
169 call convex_up_ramp_mapping_apply_backward(this%f_min, this%f_max, &
170 this%q, df_dx_in, df_dx_out, x_in)
171 else
172 call convex_down_ramp_mapping_apply_backward(this%f_min, this%f_max, &
173 this%q, df_dx_in, df_dx_out, x_in)
174 end if
175
176 end subroutine ramp_mapping_apply_backward
177
181 subroutine convex_down_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
182 real(kind=rp), intent(in) :: q, f_min, f_max
183 type(field_t), intent(in) :: x_in
184 type(field_t), intent(inout) :: X_out
185 integer :: n, i
186
187 ! x_out = f_min + (f_max - f_min) * x_in / (1 + q * (1 - x_in) )
188
189 n = x_in%dof%size()
190 if (neko_bcknd_device .eq. 1) then
191 call device_convex_down_ramp_mapping_apply(f_min, f_max, q, &
192 x_out%x_d, x_in%x_d, n)
193 else
194 do i = 1, n
195 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
196 x_in%x(i,1,1,1) / (1.0_rp + q * (1.0_rp - x_in%x(i,1,1,1) ) )
197 end do
198 end if
199
200 end subroutine convex_down_ramp_mapping_apply
201
202
207 subroutine convex_down_ramp_mapping_apply_backward(f_min, f_max, q, &
208 dF_dX_in, dF_dX_out, X_in)
209 real(kind=rp), intent(in) :: f_min, f_max, q
210 type(field_t), intent(in) :: x_in
211 type(field_t), intent(in) :: dF_dX_out
212 type(field_t), intent(inout) :: dF_dX_in
213 integer :: n, i
214
215 ! df/dx_in = df/dx_out * dx_out/dx_in
216
217 ! dx_out/dx_in = (f_min - f_max) * (q + 1) / (1 - q*(x - 1))**2
218
219 n = x_in%dof%size()
220
221 if (neko_bcknd_device .eq. 1) then
223 df_dx_in%x_d, df_dx_out%x_d, x_in%x_d, n)
224 else
225 do i = 1, n
226 df_dx_in%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) / &
227 ((1.0_rp - q * (x_in%x(i,1,1,1) - 1.0_rp))**2) * &
228 df_dx_out%x(i,1,1,1)
229 end do
230 end if
231
233
237 subroutine convex_up_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
238 real(kind=rp), intent(in) :: f_min, f_max, q
239 type(field_t), intent(in) :: x_in
240 type(field_t), intent(inout) :: X_out
241 integer :: n, i
242
243 ! x_out = f_min + (f_max - f_min) * x_in * (q + 1) / (x_in + q)
244
245 n = x_in%dof%size()
246
247 if (neko_bcknd_device .eq. 1) then
248 call device_convex_up_ramp_mapping_apply(f_min, f_max, q, &
249 x_out%x_d, x_in%x_d, n)
250 else
251 do i = 1, n
252 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
253 x_in%x(i,1,1,1) * (1.0_rp + q ) / (x_in%x(i,1,1,1) + q)
254 end do
255 end if
256
257
258 end subroutine convex_up_ramp_mapping_apply
259
260
265 subroutine convex_up_ramp_mapping_apply_backward(f_min, f_max, q, &
266 dF_dX_in, dF_dX_out, X_in)
267 real(kind=rp), intent(in) :: f_min, f_max, q
268 type(field_t), intent(in) :: x_in
269 type(field_t), intent(in) :: dF_dX_out
270 type(field_t), intent(inout) :: dF_dX_in
271 integer :: n, i
272
273 ! df/dx_in = df/dx_out * dx_out/dx_in
274
275 ! dx_out/dx_in = (f_min - f_max) * (q + 1) / (q + x)**2
276
277 n = x_in%dof%size()
278
279 if (neko_bcknd_device .eq. 1) then
281 df_dx_in%x_d, df_dx_out%x_d, x_in%x_d, n)
282 else
283 do i = 1, n
284 df_dx_in%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) / &
285 ( (x_in%x(i,1,1,1) + q)**2) * &
286 df_dx_out%x(i,1,1,1)
287 end do
288 end if
289
291end module ramp_mapping
subroutine, public device_convex_up_ramp_mapping_apply(f_min, f_max, q, x_out_d, x_in_d, n)
subroutine, public device_convex_down_ramp_mapping_apply(f_min, f_max, q, x_out_d, x_in_d, n)
subroutine, public device_convex_up_ramp_mapping_apply_backward(f_min, f_max, q, df_dx_in_d, df_dx_out_d, x_in_d, n)
subroutine, public device_convex_down_ramp_mapping_apply_backward(f_min, f_max, q, df_dx_in_d, df_dx_out_d, x_in_d, n)
Mappings to be applied to a scalar field.
Definition mapping.f90:35
A RAMP mapping of coefficients.
subroutine ramp_mapping_init_from_attributes(this, coef)
Actual constructor.
subroutine convex_up_ramp_mapping_apply(f_min, f_max, q, x_out, x_in)
Apply the mapping.
subroutine ramp_mapping_free(this)
Destructor.
subroutine convex_down_ramp_mapping_apply(f_min, f_max, q, x_out, x_in)
Apply the mapping.
subroutine ramp_mapping_init_from_json(this, json, coef)
Constructor from json.
subroutine convex_down_ramp_mapping_apply_backward(f_min, f_max, q, df_dx_in, df_dx_out, x_in)
Apply the chain rule.
subroutine ramp_mapping_apply(this, x_out, x_in)
Apply the mapping.
subroutine ramp_mapping_apply_backward(this, df_dx_in, df_dx_out, x_in)
Apply the chain rule.
subroutine convex_up_ramp_mapping_apply_backward(f_min, f_max, q, df_dx_in, df_dx_out, x_in)
Apply the chain rule.
Base abstract class for mapping.
Definition mapping.f90:44
A RAMP mapping of coefficients This is the standard RAMP described in https://doi....