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