Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
linear_mapping.f90
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 field_math, only: field_copy, field_cmult, field_cadd
38 use mapping, only: mapping_t
39 use num_types, only: rp
40 use json_module, only: json_file
41 use field, only: field_t
42 use coefs, only: coef_t
43 use json_utils, only: json_get, json_get_or_default
44 implicit none
45 private
46
48 type, public, extends(mapping_t) :: linear_mapping_t
50 real(kind=rp) :: f_min
52 real(kind=rp) :: f_max
53
54 contains
56 procedure, pass(this) :: init => linear_mapping_init_from_json
58 procedure, pass(this) :: init_from_attributes => &
59 linear_mapping_init_from_attributes
61 procedure, pass(this) :: free => linear_mapping_free
63 procedure, pass(this) :: forward_mapping => linear_forward_mapping
65 procedure, pass(this) :: backward_mapping => &
66 linear_backward_mapping
67 end type linear_mapping_t
68
69contains
70
72 subroutine linear_mapping_init_from_json(this, json, coef)
73 class(linear_mapping_t), intent(inout) :: this
74 type(json_file), intent(inout) :: json
75 type(coef_t), intent(inout) :: coef
76 real(kind=rp) :: f_min, f_max
77
78 call json_get_or_default(json, 'f_min', f_min, 0.0_rp)
79 call json_get(json, 'f_max', f_max)
80
81 call this%init_base(json, coef)
82 call linear_mapping_init_from_attributes(this, coef, f_min, f_max)
83
84 end subroutine linear_mapping_init_from_json
85
87 subroutine linear_mapping_init_from_attributes(this, coef, f_min, f_max)
88 class(linear_mapping_t), intent(inout) :: this
89 type(coef_t), intent(inout) :: coef
90 real(kind=rp), intent(in) :: f_min, f_max
91
92 this%f_min = f_min
93 this%f_max = f_max
94
95 end subroutine linear_mapping_init_from_attributes
96
98 subroutine linear_mapping_free(this)
99 class(linear_mapping_t), intent(inout) :: this
100
101 call this%free_base()
102
103 end subroutine linear_mapping_free
104
109 subroutine linear_forward_mapping(this, X_out, X_in)
110 class(linear_mapping_t), intent(inout) :: this
111 type(field_t), intent(in) :: X_in
112 type(field_t), intent(inout) :: X_out
113
114 ! x_out = f_min + (f_max - f_min) * x_in
115 call field_copy(x_out, x_in)
116 call field_cmult(x_out, this%f_max - this%f_min)
117 call field_cadd(x_out, this%f_min)
118
119 end subroutine linear_forward_mapping
120
121
127 subroutine linear_backward_mapping(this, sens_out, sens_in, X_in)
128 class(linear_mapping_t), intent(inout) :: this
129 type(field_t), intent(in) :: X_in
130 type(field_t), intent(in) :: sens_in
131 type(field_t), intent(inout) :: sens_out
132
133 ! df/dx_in = df/dx_out * dx_out/dx_in
134
135 ! dx_out/dx_in = (f_max - f_min)
136
137 call field_copy(sens_out, sens_in)
138 call field_cmult(sens_out, this%f_max - this%f_min)
139
140 end subroutine linear_backward_mapping
141
142end module linear_mapping
A linear mapping of coefficients.
Mappings to be applied to a scalar field.
Definition mapping.f90:35
A linear mapping of coefficients $f(x) = f_{min} + (f_{max} - f_{min}) x$.
Base abstract class for mapping.
Definition mapping.f90:45