Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
mask_ops.f90
Go to the documentation of this file.
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
43 use device_math, only: device_copy
44 use math_ext, only: copy_mask
45 use math, only: copy
46 use vector, only: vector_t
47 implicit none
48
49 private
51
53 module procedure mask_exterior_const_vec
54 module procedure mask_exterior_const_fld
55 end interface mask_exterior_const
56
57contains
58
63 subroutine mask_exterior_const_vec(vec, mask, const)
64 type(vector_t), intent(inout) :: vec
65 class(point_zone_t), intent(inout) :: mask
66 real(kind=rp), intent(in) :: const
67 type(field_t), pointer :: work
68 integer :: temp_indices(1), i
69
70 ! To be discussed
71 ! From what I understand with this vector/field distinction is that
72 ! ultimately the design will only contain the GLL pts inside the
73 ! the optimization domain, correct?
74 !
75 ! If this is the case, then this function makes no sense since it forces
76 ! the vector and field to be the same size.
77 !
78 ! Alternatively, this vector/field distinction is just to make the types
79 ! compatible with MMA, in which case we can continue how things are here.
80 !
81 ! In any case, it's a bit confusing and we should throw an error if the
82 ! sizes are different
83 call neko_scratch_registry%request_field(work, temp_indices(1))
84
85 if (vec%n .ne. work%size()) then
86 call neko_error('vector and field are of incompatible dimension')
87 end if
88
89 ! fill background fld
90 call field_cfill(work, const)
91
92 ! copy the fld in the masked region
93 if (neko_bcknd_device .eq. 1) then
94 call device_copy_mask(work%x_d, vec%x_d, work%size(), mask%mask_d, &
95 mask%size)
96 else
97 do i = 1, mask%size
98 work%x(mask%mask(i), 1, 1, 1) = vec%x(mask%mask(i))
99 end do
100 end if
101
102 ! copy over
103 if (neko_bcknd_device .eq. 1) then
104 call device_copy(vec%x_d, work%x_d, work%size())
105 else
106 call copy(vec%x, work%x, work%size())
107 end if
108
109 call neko_scratch_registry%relinquish_field(temp_indices)
110
111 end subroutine mask_exterior_const_vec
112
117 subroutine mask_exterior_const_fld(fld, mask, const)
118 type(field_t), intent(inout) :: fld
119 class(point_zone_t), intent(inout) :: mask
120 real(kind=rp), intent(in) :: const
121 type(field_t), pointer :: work
122 integer :: temp_indices(1)
123
124 call neko_scratch_registry%request_field(work , temp_indices(1))
125
126 ! fill background fld
127 call field_cfill(work, const)
128
129 ! copy the fld in the masked region
130 if (neko_bcknd_device .eq. 1) then
131 call device_copy_mask(work%x_d, fld%x_d, fld%size(), mask%mask_d, &
132 mask%size)
133 else
134 call copy_mask(work%x, fld%x, fld%size(), mask%mask, mask%size)
135 end if
136
137 ! copy over
138 call field_copy(fld, work)
139
140 call neko_scratch_registry%relinquish_field(temp_indices)
141
142 end subroutine mask_exterior_const_fld
143
148 subroutine mask_exterior_fld(fld, mask, background)
149 type(field_t), intent(inout) :: fld
150 class(point_zone_t), intent(inout) :: mask
151 type(field_t), intent(inout) :: background
152 type(field_t), pointer :: work
153 integer :: temp_indices(1)
154
155 call neko_scratch_registry%request_field(work , temp_indices(1))
156
157 ! fill background fld
158 call field_copy(work, background)
159
160 ! copy the fld in the masked region
161 if (neko_bcknd_device .eq. 1) then
162 call device_copy_mask(work%x_d, fld%x_d, fld%size(), mask%mask_d, &
163 mask%size)
164 else
165 call copy_mask(work%x, fld%x, fld%size(), mask%mask, mask%size)
166 end if
167
168 ! copy over
169 call field_copy(fld, work)
170
171 call neko_scratch_registry%relinquish_field(temp_indices)
172
173 end subroutine mask_exterior_fld
174end module mask_ops
subroutine device_copy_mask(a_d, b_d, size, mask_d, mask_size)
Some common Masking operations we may need.
Definition mask_ops.f90:34
subroutine, public mask_exterior_fld(fld, mask, background)
Force everything outside the mask to be a background field.
Definition mask_ops.f90:149
subroutine mask_exterior_const_vec(vec, mask, const)
Force everything outside the mask to be a constant value.
Definition mask_ops.f90:64
subroutine mask_exterior_const_fld(fld, mask, const)
Force everything outside the mask to be a constant value.
Definition mask_ops.f90:118
subroutine copy_mask(a, b, size, mask, mask_size)
Copy within a mask. NOTE, this differs from masked_copy in the indexing. .
Definition math_ext.f90:13