Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
mask_ops.f90
1! Copyright (c) 2024, 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 field, only: field_t
36 use neko_config, only: neko_bcknd_device
37 use num_types, only: rp, xp
38 use utils, only: neko_error
39 use point_zone, only: point_zone_t
40 use scratch_registry, only: neko_scratch_registry
41 use field_math, only: field_cfill, field_copy, field_rone
42 use device_math_ext, only: device_copy_mask
43 use device_math, only: device_copy, device_glsc2
44 use math_ext, only: copy_mask
45 use math, only: copy, glsc2
46 use vector, only: vector_t
47 use coefs, only: coef_t
48 implicit none
49
50 private
51 public :: mask_exterior_const, mask_exterior_fld, compute_masked_volume
52
54 module procedure mask_exterior_const_vec
55 module procedure mask_exterior_const_fld
56 end interface mask_exterior_const
57
58contains
59
64 subroutine mask_exterior_const_vec(vec, mask, const)
65 type(vector_t), intent(inout) :: vec
66 class(point_zone_t), intent(inout) :: mask
67 real(kind=rp), intent(in) :: const
68 type(field_t), pointer :: work
69 integer :: temp_indices(1), i
70
71 ! To be discussed
72 ! From what I understand with this vector/field distinction is that
73 ! ultimately the design will only contain the GLL pts inside the
74 ! the optimization domain, correct?
75 !
76 ! If this is the case, then this function makes no sense since it forces
77 ! the vector and field to be the same size.
78 !
79 ! Alternatively, this vector/field distinction is just to make the types
80 ! compatible with MMA, in which case we can continue how things are here.
81 !
82 ! In any case, it's a bit confusing and we should throw an error if the
83 ! sizes are different
84 call neko_scratch_registry%request_field(work, temp_indices(1))
85
86 if (vec%n .ne. work%size()) then
87 call neko_error('vector and field are of incompatible dimension')
88 end if
89
90 ! fill background fld
91 call field_cfill(work, const)
92
93 ! copy the fld in the masked region
94 if (neko_bcknd_device .eq. 1) then
95 call device_copy_mask(work%x_d, vec%x_d, work%size(), mask%mask_d, &
96 mask%size)
97 else
98 do i = 1, mask%size
99 work%x(mask%mask(i), 1, 1, 1) = vec%x(mask%mask(i))
100 end do
101 end if
102
103 ! copy over
104 if (neko_bcknd_device .eq. 1) then
105 call device_copy(vec%x_d, work%x_d, work%size())
106 else
107 call copy(vec%x, work%x, work%size())
108 end if
109
110 call neko_scratch_registry%relinquish_field(temp_indices)
111
112 end subroutine mask_exterior_const_vec
113
118 subroutine mask_exterior_const_fld(fld, mask, const)
119 type(field_t), intent(inout) :: fld
120 class(point_zone_t), intent(inout) :: mask
121 real(kind=rp), intent(in) :: const
122 type(field_t), pointer :: work
123 integer :: temp_indices(1)
124
125 call neko_scratch_registry%request_field(work , temp_indices(1))
126
127 ! fill background fld
128 call field_cfill(work, const)
129
130 ! copy the fld in the masked region
131 if (neko_bcknd_device .eq. 1) then
132 call device_copy_mask(work%x_d, fld%x_d, fld%size(), mask%mask_d, &
133 mask%size)
134 else
135 call copy_mask(work%x, fld%x, fld%size(), mask%mask, mask%size)
136 end if
137
138 ! copy over
139 call field_copy(fld, work)
140
141 call neko_scratch_registry%relinquish_field(temp_indices)
142
143 end subroutine mask_exterior_const_fld
144
149 subroutine mask_exterior_fld(fld, mask, background)
150 type(field_t), intent(inout) :: fld
151 class(point_zone_t), intent(inout) :: mask
152 type(field_t), intent(inout) :: background
153 type(field_t), pointer :: work
154 integer :: temp_indices(1)
155
156 call neko_scratch_registry%request_field(work , temp_indices(1))
157
158 ! fill background fld
159 call field_copy(work, background)
160
161 ! copy the fld in the masked region
162 if (neko_bcknd_device .eq. 1) then
163 call device_copy_mask(work%x_d, fld%x_d, fld%size(), mask%mask_d, &
164 mask%size)
165 else
166 call copy_mask(work%x, fld%x, fld%size(), mask%mask, mask%size)
167 end if
168
169 ! copy over
170 call field_copy(fld, work)
171
172 call neko_scratch_registry%relinquish_field(temp_indices)
173
174 end subroutine mask_exterior_fld
175
179 function compute_masked_volume(mask, coef)
180 class(point_zone_t), intent(inout) :: mask !this should be (in)
181 class(coef_t), intent(in) :: coef
182 real(kind=rp) :: compute_masked_volume
183 type(field_t), pointer :: work
184 integer :: temp_indices(1)
185 real(kind=rp) :: tmp
186 integer n
187
188 ! This would be much smarter with a kernel similar to masked_glsc2
189 ! When that kernel get written, we can update this function.
190
191 call neko_scratch_registry%request_field(work , temp_indices(1))
192 n = work%size()
193 call field_rone(work)
194 call mask_exterior_const_fld(work, mask, 0.0_rp)
195 if (neko_bcknd_device .eq. 1) then
196 tmp = device_glsc2(work%x_d, coef%B_d, n)
197 else
198 tmp = glsc2(work%x, coef%B, n)
199 end if
200 call neko_scratch_registry%relinquish_field(temp_indices)
201 compute_masked_volume = tmp
202 end function compute_masked_volume
203
204end module mask_ops
Some common Masking operations we may need.
Definition mask_ops.f90:34