Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
design_brinkman.f90
Go to the documentation of this file.
1
34
35! Implements the `brinkman_design_t` type.
36module brinkman_design
37 use num_types, only: rp, sp, dp
38 use field, only: field_t
39 use json_module, only: json_file
41 use coefs, only: coef_t
43 use scratch_registry, only: neko_scratch_registry
44 use point_zone_registry, only: neko_point_zone_registry
45 use point_zone, only: point_zone_t
47 use neko_config, only: neko_bcknd_device
48 use design, only: design_t
49 use simulation_m, only: simulation_t
51 use vector, only: vector_t
52 use matrix, only: matrix_t
53 use vector_math, only: vector_cmult
54 use matrix_math, only: matrix_cmult
55 use math, only: copy, col2
56 use device_math, only: device_copy, device_col2
57 use registry, only: neko_registry
60 use field_math, only: field_rzero
61 use json_utils, only: json_get, json_get_or_default
62 use utils, only: neko_error
63 implicit none
64 private
65
67 type, extends(design_t), public :: brinkman_design_t
68 private
69
70 ! TODO
71 ! in the future make this a derived type of a `design_variable`
72 ! type, public, extends(design_variable_t) :: brinkman_design_t
74 type(field_t), pointer :: design_indicator
75
77 ! TODO
78 ! NOTE: Tim, right now we only have the brinkman term
79 ! in the future we may map to other coeeficients for other equations...
80 ! in which case this should be a field list
81 !
82 ! or as I describe below, we also have multiple constraints,
83 ! so a list-of-lists may be the correct way forward
84 type(field_t), pointer :: brinkman_amplitude
85
86 ! NOTE:
87 ! again, we have to be so clear with nomenclature.
88 ! If we have an objective function F.
89 ! I also like to use \rho to describe the design_indicator
90 !
91 ! and let's say, for completness, we're doing conjugate heat transfer,
92 ! so we have to map
93 ! to 3 coefficients, \chi, C and \kappa.
94 !
95 ! then we will perform 3 mapping,
96 ! \rho -> \chi
97 ! \rho -> C
98 ! \rho -> \kappa
99 !
100 ! Then during the forward/adjoint looping there will be an additional
101 ! object, the "objective_function" object that will be responsible for
102 ! computing the sensitivity of the objective function with respect to the
103 ! coefficients
104 ! ie,
105 ! dF/d\chi, dF/dC and dF/d\kappa
106 !
107 ! What I'm calling "sensitivity" here, is the sensitivity with respect to
108 ! the design indicator
109 ! so dF/d\rho
110 !
111 ! so the proceedure "map_backwards" will take in the field list
112 ! dF/d\chi, dF/dC and dF/d\kappa
113 !and chain rule it's way back to
114 ! dF/d\rho
115 ! and store it here v
116 type(field_t), pointer :: sensitivity
117 ! have a field list here
118 ! type(filed_list_t), public :: constraint_sensitivity
119 ! HOWEVER !
120 ! What is a bit confusing to me... is how we'll deal with constraints.
121 !
122 ! Because in principle we could have constraints C1, C2, C3 etc
123 ! Implying we also want dC1/d\rho, dC2/d\rho etc
124 ! So this "sensitivity" should also be a list...
125 !
126 ! So I bet you'll have a nice abstract way of constructing this in the
127 ! future but I think a sort of list-of-lists would be nice.
128 !
129 ! For F:
130 ! \rho -> \tild{\rho} -> \chi
131 ! \rho -> \tild{\rho} -> C
132 ! \rho -> \tild{\rho} -> \kappa
133 !
134 ! For C1:
135 ! \rho -> \tild{\rho} (eg for a volume constraint or so)
136 !
137 ! For C2:
138 ! \rho -> \tild{\rho} (for something else)
139 !
140 ! etc..
141 !
142 ! So now we have multiple instances of the "objective" type,
143 ! each for F, C1, C2 etc
144 ! each of them can pass their dF\d_coefficents to the design
145 ! and we continue from there.
146 !
147 ! perhaps even the "objective" type is defined as a list of objectives.
148
149 ! Let's say way have a chain of two mappings
150
151 ! This is still up for discussion where the mapping will eventually live.
152 ! Maybe it makes more sense to keep the design as clean as possible,
153 ! and then the resposibility of converting the design into coefficients
154 ! used in the PDE being solved would be the responsibility of the module
155 ! used for solving the PDE, ie, the simulation.
156 ! Then the simulation would responsiblity for asking "do we have a scalar?"
157 ! ok, let me figure out how to map this design into a coeffient in the
158 ! scalar equation.
159 !
160 ! I guess then the mapping backwards gets confusing to me.
161 ! anyway, for now, let's keep things how they are and I suppose the
162 ! intention is to have different types of designs for different types of
163 ! problems.
164
170
172 class(point_zone_t), pointer :: optimization_domain
174 logical :: has_mask
176 class(coef_t), pointer :: coef
178 real(kind=rp) :: avg_b
179
180 ! TODO
181 ! you also had logicals for convergence etc,
182 ! I feel they should live in "problem"
183
184 contains
185
186 ! ----------------------------------------------------------------------- !
187 ! Initializations
188
190 generic, public :: init => init_from_json_sim, init_from_components
192 procedure, pass(this), public :: init_from_json_sim => &
193 brinkman_design_init_from_json_sim
195 procedure, pass(this), public :: init_from_components => &
196 brinkman_design_init_from_components
197
199 procedure, pass(this) :: get_values => brinkman_design_get_design
201 procedure, pass(this) :: get_sensitivity => brinkman_design_get_sensitivity
202
204 procedure, pass(this) :: design_get_x => brinkman_design_get_x
206 procedure, pass(this) :: design_get_x_i => brinkman_design_get_x_i
208 procedure, pass(this) :: design_get_y => brinkman_design_get_y
210 procedure, pass(this) :: design_get_y_i => brinkman_design_get_y_i
212 procedure, pass(this) :: design_get_z => brinkman_design_get_z
214 procedure, pass(this) :: design_get_z_i => brinkman_design_get_z_i
215
217 procedure, pass(this) :: update_design => brinkman_design_update_design
218
220 ! design_indicator -> filtering -> chi
221 ! and ultimately handle mapping different coeficients!
222 procedure, pass(this) :: map_forward => brinkman_design_map_forward
224 ! d_design_indicator <- d_filtering <- d_chi
225 ! and ultimately handle mapping different coeficients!
226 procedure, pass(this) :: map_backward => brinkman_design_map_backward
229 procedure, pass(this) :: convert_to_directional_derivative => &
230 brinkman_design_convert_to_directional_derivative
232 procedure, pass(this) :: project_sensitivity_vector => &
233 brinkman_design_project_sensitivity_vector
235 procedure, pass(this) :: project_sensitivity_matrix => &
236 brinkman_design_project_sensitivity_matrix
237 generic, public :: project_sensitivity => &
238 project_sensitivity_vector, project_sensitivity_matrix
239 ! TODO
240 ! maybe it would have been smarter to have a "coeficient" type,
241 ! which is just a scalar field and set of mappings going from
242 ! design_indicator -> coeficient and their corresponding chain rules
243 ! maybe also some information about what equation they live in...
244
245 ! a writer being called from outside would be nice
246 procedure, pass(this) :: write => brinkman_design_write
248 procedure, pass(this) :: set_output_counter => &
249 brinkman_design_set_output_counter
250
252 procedure, pass(this) :: free => brinkman_design_free
253 ! TODO
254 ! I'm not sure who owns the optimizer...
255 ! but it would make sense to have it in here so you provide it
256 ! with dF/d_design_indicator and it updates itself.
257 ! procedure, pass(this) :: update => brinkman_design_update_design
258 end type brinkman_design_t
259
260contains
261
263 subroutine brinkman_design_init_from_json_sim(this, parameters, simulation)
264 class(brinkman_design_t), intent(inout) :: this
265 type(json_file), intent(inout) :: parameters
266 type(simulation_t), intent(inout) :: simulation
267 type(json_file) :: json_subdict
268 character(len=:), allocatable :: domain_name, domain_type, name
269 character(len=:), allocatable :: output_format_str, output_precision_str
270 logical :: dealias, verbose_design, verbose_sensitivity
271 integer :: output_precision
272
273 call json_get_or_default(parameters, 'name', name, 'Brinkman Design')
274 call json_get_or_default(parameters, 'domain.type', domain_type, 'full')
275 call json_get_or_default(parameters, 'dealias', dealias, .true.)
276 call json_get_or_default(parameters, 'verbose_design', verbose_design, &
277 .false.)
278 call json_get_or_default(parameters, 'verbose_sensitivity', &
279 verbose_sensitivity, .false.)
280 call json_get_or_default(parameters, 'output_precision', &
281 output_precision_str, 'sp')
282 call json_get_or_default(parameters, 'output_format', &
283 output_format_str, 'fld')
284
285 select case (trim(output_precision_str))
286 case ('sp', 'SP', 'single')
287 output_precision = sp
288 case ('dp', 'DP', 'double')
289 output_precision = dp
290 case default
291 call neko_error('Invalid output_precision in design.brinkman. ' // &
292 'Expected "sp" or "dp".')
293 end select
294
295 select case (trim(domain_type))
296 case ('full')
297 this%has_mask = .false.
298 case ('point_zone')
299 this%has_mask = .true.
300 call json_get(parameters, 'domain.zone_name', domain_name)
301 this%optimization_domain => &
302 neko_point_zone_registry%get_point_zone(domain_name)
303
304 case default
305 call neko_error('brinkman design only supports point_zones for ' // &
306 'optimization domain types')
307
308 end select
309
310 ! Initialize and inject into the simulation
311 call this%init_from_components(name, simulation, dealias)
312
313 ! Initialize the mapper
314 associate(coef => simulation%neko_case%fluid%c_Xh, &
315 gs => simulation%neko_case%fluid%gs_Xh)
316
317 if ('mapping' .in. parameters) then
318 call this%mapping%init_base(coef)
319 call this%mapping%add(parameters, 'mapping')
320 call this%mapping%init_output_fields(this%design_indicator, &
321 this%brinkman_amplitude, this%sensitivity, verbose_design, &
322 verbose_sensitivity, output_precision, output_format_str, &
323 'design', 'sensitivity')
324 end if
325
326 if ('initial_distribution' .in. parameters) then
327 call json_get(parameters, 'initial_distribution', json_subdict)
328 call set_optimization_ic(this%design_indicator, coef, gs, &
329 json_subdict)
330 else
331 call field_rzero(this%design_indicator)
332 end if
333 end associate
334
335 ! Map to the Brinkman amplitude
336 call this%map_forward()
337
338 end subroutine brinkman_design_init_from_json_sim
339
341 subroutine brinkman_design_free(this)
342 class(brinkman_design_t), intent(inout) :: this
343
344 call this%free_base()
345 nullify(this%brinkman_amplitude)
346 nullify(this%design_indicator)
347 nullify(this%sensitivity)
348
349 end subroutine brinkman_design_free
350
351 subroutine brinkman_design_init_from_components(this, name, simulation, &
352 dealias)
353 class(brinkman_design_t), intent(inout) :: this
354 character(len=*), intent(in) :: name
355 type(simulation_t), intent(inout) :: simulation
356 logical, intent(in) :: dealias
357 integer :: n
358 type(simple_brinkman_source_term_t) :: forward_brinkman, adjoint_brinkman
359
360 associate(dof => simulation%neko_case%fluid%dm_Xh)
361
362 call neko_registry%add_field(dof, "design_indicator", .true.)
363 call neko_registry%add_field(dof, "brinkman_amplitude", .true.)
364 call neko_registry%add_field(dof, "sensitivity", .true.)
365
366 end associate
367
368 this%design_indicator => neko_registry%get_field("design_indicator")
369 this%brinkman_amplitude => neko_registry%get_field("brinkman_amplitude")
370 this%sensitivity => neko_registry%get_field("sensitivity")
371
372 this%coef => simulation%fluid%c_Xh
373
374 ! TODO
375 ! this is where we steal basically everything in
376 ! brinkman_source_term regarding loading initial fields
377 ! for now, make it a cylinder by hand
378 this%design_indicator = 0.0_rp
379 this%brinkman_amplitude = 0.0_rp
380 this%design_indicator = 0.0_rp
381
382 ! TODO
383 ! Regarding masks and filters,
384 ! I suppose there are two ways of thinking about it:
385 ! 1) Mask first, then filter
386 ! 2) filter first, then mask
387 !
388 ! Each one can be used to achieve different results, and when we do complex
389 ! non-linear filter cascades the choice here has implications for exactly
390 ! how we define "minimum size control".
391 !
392 ! On top of this, I've only ever used spatial convolution filters before,
393 ! not PDE based filters.
394 ! With these filters you need to think of how your domain is "padded" when
395 ! you move the kernel over the boundary. Again, you can achieve different
396 ! effects depending on the choice of padding.
397 ! I don't really know if there's a similar implication for PDE based
398 ! based filters, I bet there is. (We should ask Niels and Casper) I bet it
399 ! means we need to consider the boundary of the mask as the boundary we
400 ! enforce to be Nuemann, not the boundary of the computational domain.
401 ! That sounds REALLY hard to implement...I hope it doesn't come down to
402 ! that.
403 !
404 ! We can always change this decision later, but I'm going with (1)
405 ! mask first, then filter.
406 ! The reason being, one of the purposes of the filtering is to avoid sharp
407 ! interfaces, if we filter first then mask there's a chance we have a sharp
408 ! interface on the boundary of the optimization domain.
409 ! if we mask first then filter, at least all the boundaries will be smooth.
410 if (this%has_mask) then
411 call mask_exterior_const(this%design_indicator, &
412 this%optimization_domain, 0.0_rp)
413 end if
414
415 n = this%design_indicator%dof%size()
416 call this%init_base(name, n)
417
418 ! compute the average mass matrix
419 this%avg_B = this%coef%volume / real(simulation%fluid%glb_unique_points)
420
421 ! init the simple brinkman term for the forward problem
422 call forward_brinkman%init_from_components( &
423 simulation%fluid%f_x, &
424 simulation%fluid%f_y, &
425 simulation%fluid%f_z, &
426 this%brinkman_amplitude, &
427 simulation%fluid%u, &
428 simulation%fluid%v, &
429 simulation%fluid%w, &
430 simulation%fluid%c_Xh, &
431 simulation%adjoint_fluid%c_Xh_GL, &
432 simulation%adjoint_fluid%GLL_to_GL, &
433 dealias, simulation%adjoint_fluid%scratch_GL)
434 ! append brinkman source term to the forward problem
435 call simulation%fluid%source_term%add(forward_brinkman)
436
437 ! init the simple brinkman term for the adjoint
438 call adjoint_brinkman%init_from_components( &
439 simulation%adjoint_fluid%f_adj_x, &
440 simulation%adjoint_fluid%f_adj_y, &
441 simulation%adjoint_fluid%f_adj_z, &
442 this%brinkman_amplitude, &
443 simulation%adjoint_fluid%u_adj, &
444 simulation%adjoint_fluid%v_adj, &
445 simulation%adjoint_fluid%w_adj, &
446 simulation%adjoint_fluid%c_Xh, &
447 simulation%adjoint_fluid%c_Xh_GL, &
448 simulation%adjoint_fluid%GLL_to_GL, &
449 dealias, simulation%adjoint_fluid%scratch_GL)
450 ! append brinkman source term based on design
451
452 select type (f => simulation%adjoint_fluid)
453 type is (adjoint_fluid_pnpn_t)
454 call f%source_term%add(adjoint_brinkman)
455 class default
456 end select
457
458 end subroutine brinkman_design_init_from_components
459
460
461 subroutine brinkman_design_map_forward(this)
462 class(brinkman_design_t), intent(inout) :: this
463
464 ! TODO, see previous todo about mask first, then mapping
465 if (this%has_mask) then
466 call mask_exterior_const(this%design_indicator, &
467 this%optimization_domain, 0.0_rp)
468 end if
469
470 call this%mapping%apply_forward(this%brinkman_amplitude, &
471 this%design_indicator)
472
473 end subroutine brinkman_design_map_forward
474
475 subroutine brinkman_design_get_design(this, values)
476 class(brinkman_design_t), intent(in) :: this
477 type(vector_t), intent(inout) :: values
478 integer :: n
479
480 n = this%size()
481 if (n .ne. values%size()) then
482 call neko_error('Get design: size mismatch')
483 end if
484
485 if (neko_bcknd_device .eq. 1) then
486 call device_copy(values%x_d, this%design_indicator%x_d, n)
487 else
488 call copy(values%x, this%design_indicator%x, n)
489 end if
490
491 end subroutine brinkman_design_get_design
492
493 subroutine brinkman_design_get_sensitivity(this, values)
494 class(brinkman_design_t), intent(in) :: this
495 type(vector_t), intent(inout) :: values
496 integer :: n
497
498 n = this%size()
499 if (n .ne. values%size()) then
500 call neko_error('Get sensitivity: size mismatch')
501 end if
502
503 if (neko_bcknd_device .eq. 1) then
504 call device_copy(values%x_d, this%sensitivity%x_d, n)
505 else
506 call copy(values%x, this%sensitivity%x, n)
507 end if
508
509 end subroutine brinkman_design_get_sensitivity
510
511 ! The interpretation of the sensitivity DF is a function that satisfies
512 ! <DF,phi> = delta F(x;phi)
513 ! To recover directional derivatives we should project the sensitivity onto
514 ! the basis functions, equivelent to multiplying by the mass matrix.
515 subroutine brinkman_design_convert_to_directional_derivative(this, vec_in)
516 class(brinkman_design_t), intent(in) :: this
517 type(vector_t), intent(inout) :: vec_in
518 integer :: n
519
520 n = this%size()
521 if (n .ne. vec_in%size()) then
522 call neko_error('convert_to_directional_derivative: size mismatch')
523 end if
524
525 if (neko_bcknd_device .eq. 1) then
526 call device_col2(vec_in%x_d, this%coef%B_d, n)
527 else
528 call col2(vec_in%x, this%coef%B, n)
529 end if
530
531 end subroutine brinkman_design_convert_to_directional_derivative
532
533 ! to aleviate lumpy designs due to non uniform spacing, we instead project
534 ! onto a uniform mass matrix with average weight, as opposed to using the
535 ! directional derivative
536 subroutine brinkman_design_project_sensitivity_vector(this, sensitivity)
537 class(brinkman_design_t), intent(inout) :: this
538 type(vector_t), intent(inout) :: sensitivity
539
540 call vector_cmult(sensitivity, this%avg_B)
541 end subroutine brinkman_design_project_sensitivity_vector
542
543 subroutine brinkman_design_project_sensitivity_matrix(this, sensitivity)
544 class(brinkman_design_t), intent(inout) :: this
545 type(matrix_t), intent(inout) :: sensitivity
546
547 call matrix_cmult(sensitivity, this%avg_B)
548 end subroutine brinkman_design_project_sensitivity_matrix
549
550 subroutine brinkman_design_get_x(this, x)
551 class(brinkman_design_t), intent(in) :: this
552 type(vector_t), intent(inout) :: x
553 integer :: n
554
555 n = this%size()
556 if (n .ne. x%size()) then
557 call neko_error('Get x: size mismatch')
558 end if
559
560 if (neko_bcknd_device .eq. 1) then
561 call device_copy(x%x_d, this%design_indicator%dof%x_d, n)
562 else
563 call copy(x%x, this%design_indicator%dof%x, n)
564 end if
565
566 end subroutine brinkman_design_get_x
567
568 function brinkman_design_get_x_i(this, i) result(x_i)
569 class(brinkman_design_t), intent(in) :: this
570 integer, intent(in) :: i
571 real(kind=rp) :: x_i
572 integer :: n
573
574 n = this%size()
575 if (i .lt. 1 .or. i .gt. n) then
576 call neko_error('brinkman_design_get_x_i: index out of bounds')
577 end if
578
579 x_i = this%design_indicator%dof%x(i,1,1,1)
580
581 end function brinkman_design_get_x_i
582
583 subroutine brinkman_design_get_y(this, y)
584 class(brinkman_design_t), intent(in) :: this
585 type(vector_t), intent(inout) :: y
586 integer :: n
587
588 n = this%size()
589 if (n .ne. y%size()) then
590 call neko_error('Get y: size mismatch')
591 end if
592
593 if (neko_bcknd_device .eq. 1) then
594 call device_copy(y%x_d, this%design_indicator%dof%y_d, n)
595 else
596 call copy(y%x, this%design_indicator%dof%y, n)
597 end if
598
599 end subroutine brinkman_design_get_y
600
601 function brinkman_design_get_y_i(this, i) result(y_i)
602 class(brinkman_design_t), intent(in) :: this
603 integer, intent(in) :: i
604 real(kind=rp) :: y_i
605 integer :: n
606
607 n = this%size()
608 if (i .lt. 1 .or. i .gt. n) then
609 call neko_error('brinkman_design_get_y_i: index out of bounds')
610 end if
611
612 y_i = this%design_indicator%dof%y(i,1,1,1)
613
614 end function brinkman_design_get_y_i
615
616 subroutine brinkman_design_get_z(this, z)
617 class(brinkman_design_t), intent(in) :: this
618 type(vector_t), intent(inout) :: z
619 integer :: n
620
621 n = this%size()
622 if (n .ne. z%size()) then
623 call neko_error('Get z: size mismatch')
624 end if
625
626 if (neko_bcknd_device .eq. 1) then
627 call device_copy(z%x_d, this%design_indicator%dof%z_d, n)
628 else
629 call copy(z%x, this%design_indicator%dof%z, n)
630 end if
631
632 end subroutine brinkman_design_get_z
633
634 function brinkman_design_get_z_i(this, i) result(z_i)
635 class(brinkman_design_t), intent(in) :: this
636 integer, intent(in) :: i
637 real(kind=rp) :: z_i
638 integer :: n
639
640 n = this%size()
641 if (i .lt. 1 .or. i .gt. n) then
642 call neko_error('brinkman_design_get_z_i: index out of bounds')
643 end if
644
645 z_i = this%design_indicator%dof%z(i,1,1,1)
646
647 end function brinkman_design_get_z_i
648
649 subroutine brinkman_design_update_design(this, values)
650 class(brinkman_design_t), intent(inout) :: this
651 type(vector_t), intent(inout) :: values
652 integer :: n
653
654 n = this%size()
655 if (neko_bcknd_device .eq. 1) then
656 call device_copy(this%design_indicator%x_d, values%x_d, n)
657 else
658 call copy(this%design_indicator%x, values%x, n)
659 end if
660
661 call this%map_forward()
662
663 if (neko_bcknd_device .eq. 1) then
664 call device_copy(values%x_d, this%design_indicator%x_d, n)
665 else
666 call copy(values%x, this%design_indicator%x, n)
667 end if
668
669 end subroutine brinkman_design_update_design
670
671 subroutine brinkman_design_map_backward(this, sensitivity)
672 class(brinkman_design_t), intent(inout) :: this
673 type(vector_t), intent(in) :: sensitivity
674 type(field_t), pointer :: tmp_fld
675 integer :: temp_index
676
677 call neko_scratch_registry%request(tmp_fld, temp_index, .false.)
678
679 call vector_to_field(tmp_fld, sensitivity)
680
681 call this%mapping%apply_backward(this%sensitivity, tmp_fld)
682
683 if (this%has_mask) then
684 call mask_exterior_const(this%sensitivity, this%optimization_domain, &
685 0.0_rp)
686 end if
687
688 call neko_scratch_registry%relinquish_field(temp_index)
689
690 end subroutine brinkman_design_map_backward
691
692 subroutine brinkman_design_write(this, idx)
693 class(brinkman_design_t), intent(inout) :: this
694 integer, intent(in) :: idx
695
696 call this%mapping%write_design(idx)
697 call this%mapping%write_sensitivity(idx)
698
699 end subroutine brinkman_design_write
700
701 subroutine brinkman_design_set_output_counter(this, idx)
702 class(brinkman_design_t), intent(inout) :: this
703 integer, intent(in) :: idx
704
705 call this%mapping%set_output_counter(idx)
706
707 end subroutine brinkman_design_set_output_counter
708
709end module brinkman_design
Adjoint Pn/Pn formulation.
Implements the design_t.
Definition design.f90:36
Implements the mapping_handler_t type.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
Contains extensions to the neko library required to run the topology optimization code.
Definition neko_ext.f90:42
subroutine, public field_to_vector(vector, field)
Field to vector.
Definition neko_ext.f90:409
subroutine, public vector_to_field(field, vector)
Vector to field.
Definition neko_ext.f90:385
Optimization initial condition.
Implements the simple_brinkman_source_term_t type.
Implements the steady_problem_t type.
A topology optimization design variable.
An abstract design type.
Definition design.f90:53
Abstract class for handling mapping_cascade.