V3FIT
limiter.f
Go to the documentation of this file.
1 !*******************************************************************************
4 !
5 ! Note separating the Doxygen comment block here so detailed decription is
6 ! found in the Module not the file.
7 !
13 !*******************************************************************************
14 
15  MODULE limiter
16 
17  USE stel_kinds, only: rprec, dp
18  USE profiler
19  USE signal
20 
21  IMPLICIT NONE
22 !*******************************************************************************
23 ! geometric module parameters
24 !*******************************************************************************
25 !*******************************************************************************
26 ! DERIVED-TYPE DECLARATIONS
27 ! 1) limiter base class
28 !
29 !*******************************************************************************
30 !-------------------------------------------------------------------------------
34 !-------------------------------------------------------------------------------
35  TYPE, EXTENDS(signal_class) :: limiter_class
38  LOGICAL :: on_edge = .false.
40  REAL (rprec), DIMENSION(:), POINTER :: phi => null()
41  CONTAINS
42  PROCEDURE :: &
43  & get_modeled_signal_last => limiter_get_modeled_signal
44  PROCEDURE ::
45  & get_header => limiter_get_header
46  PROCEDURE :: &
47  & get_type => limiter_get_type
48  PROCEDURE :: &
49  & get_max_fval => limiter_get_max_fval
50  END TYPE
51 
52 !*******************************************************************************
53 ! INTERFACE BLOCKS
54 !*******************************************************************************
55 
56  CONTAINS
57 
58 !*******************************************************************************
59 ! DESTRUCTION SUBROUTINES
60 !*******************************************************************************
61 !-------------------------------------------------------------------------------
67 !-------------------------------------------------------------------------------
68  SUBROUTINE limiter_destruct(this)
69 
70  IMPLICIT NONE
71 
72 ! Declare Arguments
73  class(limiter_class), INTENT(inout) :: this
74 
75 ! Start of executable code
76  IF (ASSOCIATED(this%phi)) THEN
77  DEALLOCATE(this%phi)
78  this%phi => null()
79  END IF
80 
81  this%on_edge = .false.
82 
83  END SUBROUTINE
84 
85 !*******************************************************************************
86 ! GETTER SUBROUTINES
87 !*******************************************************************************
88 !-------------------------------------------------------------------------------
102 !-------------------------------------------------------------------------------
103  FUNCTION limiter_get_modeled_signal(this, a_model, sigma, &
104  & last_value)
105  USE model
106 
107  IMPLICIT NONE
108 
109 ! Declare Arguments
110  REAL (rprec), DIMENSION(4) :: limiter_get_modeled_signal
111  class(limiter_class), INTENT(inout) :: this
112  TYPE (model_class), POINTER :: a_model
113  REAL (rprec), DIMENSION(4), INTENT(out) :: sigma
114  REAL (rprec), DIMENSION(4), INTENT(in) :: last_value
115 
116 ! local variables
117 ! The r and z arrays will be allocated by the equilibrium.
118  REAL (rprec), DIMENSION(:), POINTER :: r
119  REAL (rprec), DIMENSION(:), POINTER :: z
120  INTEGER :: num_theta
121  INTEGER :: phi_index
122  REAL (rprec), DIMENSION(3) :: rphiz_at_max
123  REAL (rprec) :: fval_max
124  REAL (rprec) :: start_time
125 
126 ! Start of executable code
127  start_time = profiler_get_start_time()
128 
129  sigma = 0.0
130 
131  IF (btest(a_model%state_flags, model_state_vmec_flag) .or. &
132  & btest(a_model%state_flags, model_state_siesta_flag) .or. &
133  & btest(a_model%state_flags, model_state_shift_flag) .or. &
134  & btest(a_model%state_flags, model_state_signal_flag)) THEN
135  r => null()
136  z => null()
137 
138 ! Loop over each phi plane. r, z are allocated by equilibrium_get_plasma_edge.
139 ! There is no need to allocate and deallocate r and z for each phi angle.
140  DO phi_index = 1, SIZE(this%phi)
141  num_theta = &
142  & equilibrium_get_plasma_edge(a_model%equilibrium, &
143  & this%phi(phi_index), r, z)
144 
145  fval_max = this%get_max_fval(num_theta, phi_index, r, z, &
146  & rphiz_at_max)
147  END DO
148 
149  IF (this%on_edge) THEN
150  limiter_get_modeled_signal(1) = fval_max
151  ELSE
152  limiter_get_modeled_signal(1) = max(fval_max, 0.0)
153  END IF
154  limiter_get_modeled_signal(2:4) = rphiz_at_max
155 
156 ! Dealloctae r and z since they are no longer needed.
157  IF (ASSOCIATED(r)) THEN
158  DEALLOCATE(r)
159  END IF
160  IF (ASSOCIATED(z)) THEN
161  DEALLOCATE(z)
162  END IF
163  r => null()
164  z => null()
165 
166  CALL this%scale_and_offset(a_model, &
168  ELSE
169  limiter_get_modeled_signal = last_value
170  END IF
171 
172  CALL profiler_set_stop_time('limiter_get_modeled_signal', &
173  & start_time)
174 
175  END FUNCTION
176 
177 !-------------------------------------------------------------------------------
184 !-------------------------------------------------------------------------------
185  FUNCTION limiter_get_type(this)
187 
188  IMPLICIT NONE
189 
190 ! Declare Arguments
191  CHARACTER (len=data_name_length) :: limiter_get_type
192  class(limiter_class), INTENT(in) :: this
193 
194 ! local variables
195  REAL (rprec) :: start_time
196 
197 ! Start of executable code
198  start_time = profiler_get_start_time()
199 
200  limiter_get_type = 'edge_lim '
201 
202  CALL profiler_set_stop_time('limiter_get_type', start_time)
203 
204  END FUNCTION
205 
206 !-------------------------------------------------------------------------------
214 !-------------------------------------------------------------------------------
215  SUBROUTINE limiter_get_header(this, header)
217 
218  IMPLICIT NONE
219 
220 ! Declare Arguments
221  class(limiter_class), INTENT(in) :: this
222  CHARACTER (len=data_name_length), DIMENSION(7), INTENT(inout) :: &
223  & header
224 
225 ! local variables
226  REAL (rprec) :: start_time
227 
228 ! Start of executable code
229  start_time = profiler_get_start_time()
230 
231  header(1) = 'r (m)'
232  header(2) = 'phi (rad)'
233  header(3) = 'z (m)'
234  header(4) = 'model_sig(1)'
235  header(5) = 'model_sig(2)'
236  header(6) = 'model_sig(3)'
237  header(7) = 'model_sig(4)'
238 
239  CALL profiler_set_stop_time('limiter_get_header', start_time)
240 
241  END SUBROUTINE
242 
243 !-------------------------------------------------------------------------------
255 !-------------------------------------------------------------------------------
256  FUNCTION limiter_get_max_fval(this, num_theta, phi_index, &
257  & r, z, rphiz_at_max)
258  USE v3_utilities
259 
260  IMPLICIT NONE
261 
262 ! Declare Arguments
263  REAL (rprec) :: limiter_get_max_fval
264  CLASS (limiter_class), INTENT(in) :: this
265  INTEGER, INTENT(in) :: num_theta
266  INTEGER, INTENT(in) :: phi_index
267  REAL (rprec), DIMENSION(:), INTENT(in) :: r
268  REAL (rprec), DIMENSION(:), INTENT(in) :: z
269  REAL (rprec), DIMENSION(3), INTENT(out) :: rphiz_at_max
270 
271 ! Start of executable code
272  CALL assert(.false., 'limiter_get_max_fval not over written' // &
273  & ' for ' // this%get_type())
274 
275  rphiz_at_max = 0.0
277 
278  END FUNCTION
279 
280  END MODULE
profiler
Defines functions for measuring an tabulating performance of function and subroutine calls....
Definition: profiler.f:13
limiter::limiter_get_header
subroutine limiter_get_header(this, header)
Gets a discription of the model and model sigma array indices.
Definition: limiter.f:216
limiter::limiter_class
Base class representing a limiter signal.
Definition: limiter.f:35
model
Defines the base class of the type model_class. The model contains information not specific to the eq...
Definition: model.f:59
limiter::limiter_destruct
subroutine limiter_destruct(this)
Deconstruct a limiter_class object.
Definition: limiter.f:69
v3_utilities::assert
Definition: v3_utilities.f:55
limiter::limiter_get_modeled_signal
real(rprec) function, dimension(4) limiter_get_modeled_signal(this, a_model, sigma, last_value)
Calculates the modeled signal.
Definition: limiter.f:105
profiler::profiler_get_start_time
real(rprec) function profiler_get_start_time()
Gets the start time of profiled function.
Definition: profiler.f:194
model::model_class
Base class representing a model.
Definition: model.f:141
limiter::limiter_get_max_fval
real(rprec) function limiter_get_max_fval(this, num_theta, phi_index, r, z, rphiz_at_max)
Calculates the maximum value of the limiter function.
Definition: limiter.f:258
limiter::limiter_get_type
character(len=data_name_length) function limiter_get_type(this)
Gets a discription of the limiter type.
Definition: limiter.f:186
data_parameters
This modules contains parameters used by equilibrium models.
Definition: data_parameters.f:10
profiler::profiler_set_stop_time
subroutine profiler_set_stop_time(symbol_name, start_time)
Gets the end time of profiled function.
Definition: profiler.f:121
limiter
Defines the base class of the type limiter_class.
Definition: limiter.f:15
signal::signal_class
Base class representing a signal.
Definition: signal.f:33
signal
Defines the base class of the type signal_class.
Definition: signal.f:14