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