36 use num_types,
only : rp
37 use field_list,
only : field_list_t
38 use field,
only: field_t
39 use field_registry,
only: neko_field_registry
40 use scratch_registry,
only: neko_scratch_registry
41 use json_module,
only : json_file
42 use json_utils,
only: json_get, json_get_or_default
43 use source_term,
only : source_term_t
44 use coefs,
only : coef_t
45 use field_math,
only: field_add2s2, field_copy, field_cadd
47 use point_zone,
only: point_zone_t
62 type(field_t),
pointer :: s => null()
64 real(kind=rp) :: obj_scale
66 real(kind=rp) :: phi_ref
68 class(point_zone_t),
pointer :: mask => null()
72 real(kind=rp) :: mask_volume
75 procedure, pass(this) :: init => &
76 adjoint_mixing_scalar_source_term_init_from_json
78 procedure, pass(this) :: init_from_components => &
79 adjoint_mixing_scalar_source_term_init_from_components
81 procedure, pass(this) :: free => adjoint_mixing_scalar_source_term_free
83 procedure, pass(this) :: compute_ => &
84 adjoint_mixing_scalar_source_term_compute
94 subroutine adjoint_mixing_scalar_source_term_init_from_json(this, &
95 json, fields, coef, variable_name)
97 type(json_file),
intent(inout) :: json
98 type(field_list_t),
intent(in),
target :: fields
99 type(coef_t),
intent(in),
target :: coef
100 character(len=*),
intent(in) :: variable_name
103 end subroutine adjoint_mixing_scalar_source_term_init_from_json
106 subroutine adjoint_mixing_scalar_source_term_init_from_components(this, &
107 f_s, s, obj_scale, phi_ref, mask, if_mask, coef)
109 type(field_t),
pointer,
intent(in) :: f_s
110 real(kind=rp),
intent(in) :: phi_ref
111 type(field_t),
intent(in),
target :: s
112 class(point_zone_t),
intent(in),
target :: mask
113 type(field_list_t) :: fields
115 real(kind=rp) :: start_time
116 real(kind=rp) :: end_time
117 real(kind=rp) :: obj_scale
122 end_time = 100000000.0_rp
130 call fields%assign(1, f_s)
132 call this%init_base(fields, coef, start_time, end_time)
136 this%obj_scale = obj_scale
137 this%phi_ref = phi_ref
138 this%if_mask = if_mask
140 if (this%if_mask)
then
142 this%mask_volume = compute_masked_volume(this%mask, coef)
144 this%mask_volume = coef%volume
147 end subroutine adjoint_mixing_scalar_source_term_init_from_components
150 subroutine adjoint_mixing_scalar_source_term_free(this)
153 call this%free_base()
154 end subroutine adjoint_mixing_scalar_source_term_free
160 subroutine adjoint_mixing_scalar_source_term_compute(this, t, tstep)
162 real(kind=rp),
intent(in) :: t
163 integer,
intent(in) :: tstep
164 type(field_t),
pointer :: fs
165 type(field_t),
pointer :: work
166 integer :: temp_indices(1)
169 fs => this%fields%get(1)
171 call neko_scratch_registry%request_field(work, temp_indices(1))
173 call field_copy(work, this%s)
175 call field_cadd(work, -this%phi_ref)
177 if (this%if_mask)
then
182 call field_add2s2(fs, work, this%obj_scale / this%mask_volume)
183 call neko_scratch_registry%relinquish_field(temp_indices)
186 end subroutine adjoint_mixing_scalar_source_term_compute
Implements the adjoint_mixing_scalar_source_term type.
Some common Masking operations we may need.