source: palm/trunk/SOURCE/wind_turbine_model_mod.f90 @ 2322

Last change on this file since 2322 was 2322, checked in by Giersch, 7 years ago

Bugfix of error message and assign error numbers

  • Property svn:keywords set to Id
File size: 109.1 KB
RevLine 
[1820]1!> @file wind_turbine_model_mod.f90
[1912]2!------------------------------------------------------------------------------!
[1819]3! This file is part of PALM.
4!
[2000]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
[1819]9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2257]17! Copyright 2009-2017 Carl von Ossietzky Universitaet Oldenburg
[2101]18! Copyright 1997-2017 Leibniz Universitaet Hannover
[1912]19!------------------------------------------------------------------------------!
[1819]20!
21! Current revisions:
22! -----------------
[1913]23!
[2233]24!
[1913]25! Former revisions:
26! -----------------
27! $Id: wind_turbine_model_mod.f90 2322 2017-07-26 08:30:28Z Giersch $
[2322]28! Bugfix of error message and assign error numbers
29!
30! 2257 2017-06-07 14:07:05Z witha
[2257]31! Bugfix: turb_cl_tab and turb_cd_tab were set to zero before being allocated
[1844]32!
[2257]33! 2233 2017-05-30 18:08:54Z suehring
34!
[2233]35! 2232 2017-05-30 17:47:52Z suehring
36! Adjustments to new topography concept
37!
[2153]38! 2152 2017-02-17 13:27:24Z lvollmer
39! Bugfix in subroutine wtm_read_blade_tables
40! Addition of a tip loss model
41!
[2016]42! 2015 2016-09-28 08:45:18Z lvollmer
43! Bugfix of pitch control
44!
[2001]45! 2000 2016-08-20 18:09:15Z knoop
46! Forced header and separation lines into 80 columns
47!
[1930]48! 1929 2016-06-09 16:25:25Z suehring
49! Bugfix: added preprocessor directives for parallel and serial mode
50!
[1917]51! 1914 2016-05-26 14:44:07Z witha
[1821]52! Initial revision
53!
[1819]54!
55! Description:
56! ------------
57!> This module calculates the effect of wind turbines on the flow fields. The
58!> initial version contains only the advanced actuator disk with rotation method
59!> (ADM-R).
60!> The wind turbines include the tower effect, can be yawed and tilted.
61!> The wind turbine model includes controllers for rotational speed, pitch and
62!> yaw.
63!> Currently some specifications of the NREL 5 MW reference turbine
64!> are hardcoded whereas most input data comes from separate files (currently
65!> external, planned to be included as namelist which will be read in
66!> automatically).
67!>
68!> @todo Revise code according to PALM Coding Standard
69!> @todo Implement ADM and ALM turbine models
70!> @todo Generate header information
[1917]71!> @todo Implement further parameter checks and error messages
[1819]72!> @todo Revise and add code documentation
73!> @todo Output turbine parameters as timeseries
74!> @todo Include additional output variables
[1864]75!> @todo Revise smearing the forces for turbines in yaw
76!> @todo Revise nacelle and tower parameterization
77!> @todo Allow different turbine types in one simulation
[1819]78!
79!------------------------------------------------------------------------------!
80 MODULE wind_turbine_model_mod
81
82    USE arrays_3d,                                                             &
83        ONLY:  dd2zu, tend, u, v, w, zu, zw
84
85    USE constants
86
87    USE control_parameters,                                                    &
88        ONLY:  dt_3d, dz, message_string, simulated_time
89
90    USE cpulog,                                                                &
91        ONLY:  cpu_log, log_point_s
92
93    USE grid_variables,                                                        &
94        ONLY:  ddx, dx, ddy, dy
95
96    USE indices,                                                               &
97        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nz,   &
[2232]98               nzb, nzt, wall_flags_0
[1819]99
100    USE kinds
101
102    USE pegrid
103
104
105    IMPLICIT NONE
106
[1864]107    PRIVATE
[1819]108
[1839]109    LOGICAL ::  wind_turbine=.FALSE.   !< switch for use of wind turbine model
[1819]110
111!
112!-- Variables specified in the namelist wind_turbine_par
113
[1864]114    INTEGER(iwp) ::  nairfoils = 8   !< number of airfoils of the used turbine model (for ADM-R and ALM)
[1839]115    INTEGER(iwp) ::  nturbines = 1   !< number of turbines
[1819]116
[1839]117    LOGICAL ::  pitch_control = .FALSE.   !< switch for use of pitch controller
118    LOGICAL ::  speed_control = .FALSE.   !< switch for use of speed controller
119    LOGICAL ::  yaw_control   = .FALSE.   !< switch for use of yaw controller
[2152]120    LOGICAL ::  tl_cor        = .FALSE.    !< switch for use of tip loss correct.
[1819]121
[1839]122    REAL(wp) ::  segment_length  = 1.0_wp          !< length of the segments, the rotor area is divided into
123                                                   !< (in tangential direction, as factor of MIN(dx,dy,dz))
124    REAL(wp) ::  segment_width   = 0.5_wp          !< width of the segments, the rotor area is divided into
125                                                   !< (in radial direction, as factor of MIN(dx,dy,dz))
126    REAL(wp) ::  time_turbine_on = 0.0_wp          !< time at which turbines are started
127    REAL(wp) ::  tilt            = 0.0_wp          !< vertical tilt of the rotor [degree] ( positive = backwards )
[1819]128
[1912]129    REAL(wp), DIMENSION(1:100) ::  dtow             = 0.0_wp  !< tower diameter [m]
130    REAL(wp), DIMENSION(1:100) ::  omega_rot        = 0.0_wp  !< inital or constant rotor speed [rad/s]
131    REAL(wp), DIMENSION(1:100) ::  phi_yaw          = 0.0_wp  !< yaw angle [degree] ( clockwise, 0 = facing west )
132    REAL(wp), DIMENSION(1:100) ::  pitch_add        = 0.0_wp  !< constant pitch angle
133    REAL(wp), DIMENSION(1:100) ::  rcx        = 9999999.9_wp  !< position of hub in x-direction
134    REAL(wp), DIMENSION(1:100) ::  rcy        = 9999999.9_wp  !< position of hub in y-direction
135    REAL(wp), DIMENSION(1:100) ::  rcz        = 9999999.9_wp  !< position of hub in z-direction
136    REAL(wp), DIMENSION(1:100) ::  rnac             = 0.0_wp  !< nacelle diameter [m]
137    REAL(wp), DIMENSION(1:100) ::  rr              = 63.0_wp  !< rotor radius [m]
138    REAL(wp), DIMENSION(1:100) ::  turb_cd_nacelle = 0.85_wp  !< drag coefficient for nacelle
139    REAL(wp), DIMENSION(1:100) ::  turb_cd_tower    = 1.2_wp  !< drag coefficient for tower
[1839]140
[1819]141!
142!-- Variables specified in the namelist for speed controller
143!-- Default values are from the NREL 5MW research turbine (Jonkman, 2008)
144
[1912]145    REAL(wp) ::  rated_power    = 5296610.0_wp    !< rated turbine power [W]
[1839]146    REAL(wp) ::  gear_ratio     = 97.0_wp         !< Gear ratio from rotor to generator
147    REAL(wp) ::  inertia_rot    = 34784179.0_wp   !< Inertia of the rotor [kg/m2]
148    REAL(wp) ::  inertia_gen    = 534.116_wp      !< Inertia of the generator [kg/m2]
149    REAL(wp) ::  gen_eff        = 0.944_wp        !< Electric efficiency of the generator
150    REAL(wp) ::  gear_eff       = 1.0_wp          !< Loss between rotor and generator
151    REAL(wp) ::  air_dens       = 1.225_wp        !< Air density to convert to W [kg/m3]
152    REAL(wp) ::  rated_genspeed = 121.6805_wp     !< Rated generator speed [rad/s]
[1912]153    REAL(wp) ::  max_torque_gen = 47402.91_wp     !< Maximum of the generator torque [Nm]
[1839]154    REAL(wp) ::  slope2         = 2.332287_wp     !< Slope constant for region 2
155    REAL(wp) ::  min_reg2       = 91.21091_wp     !< Lower generator speed boundary of region 2 [rad/s]
156    REAL(wp) ::  min_reg15      = 70.16224_wp     !< Lower generator speed boundary of region 1.5 [rad/s]
[1912]157    REAL(wp) ::  max_trq_rate   = 15000.0_wp      !< Max generator torque increase [Nm/s]
158    REAL(wp) ::  pitch_rate     = 8.0_wp          !< Max pitch rate [degree/s]
[1839]159
[1912]160
[1819]161!
162!-- Variables specified in the namelist for yaw control
163
[1839]164    REAL(wp) ::  yaw_speed = 0.005236_wp   !< speed of the yaw actuator [rad/s]
165    REAL(wp) ::  max_miss = 0.08726_wp     !< maximum tolerated yaw missalignment [rad]
166    REAL(wp) ::  min_miss = 0.008726_wp    !< minimum yaw missalignment for which the actuator stops [rad]
167
[1819]168!
169!-- Set flag for output files TURBINE_PARAMETERS
170    TYPE file_status
171       LOGICAL ::  opened, opened_before
172    END TYPE file_status
173   
[1912]174    TYPE(file_status), DIMENSION(500) :: openfile_turb_mod =                   &
175                                         file_status(.FALSE.,.FALSE.)
[1819]176
177!
178!-- Variables for initialization of the turbine model
179
[1839]180    INTEGER(iwp) ::  inot         !< turbine loop index (turbine id)
181    INTEGER(iwp) ::  nsegs_max    !< maximum number of segments (all turbines, required for allocation of arrays)
182    INTEGER(iwp) ::  nrings_max   !< maximum number of rings (all turbines, required for allocation of arrays)
183    INTEGER(iwp) ::  ring         !< ring loop index (ring number)
184    INTEGER(iwp) ::  rr_int       !<
185    INTEGER(iwp) ::  upper_end    !<
[1819]186
[1839]187    INTEGER(iwp), DIMENSION(1) ::  lct   !<
[1819]188
[1912]189    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  i_hub     !< index belonging to x-position of the turbine
190    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  i_smear   !< index defining the area for the smearing of the forces (x-direction)
191    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  j_hub     !< index belonging to y-position of the turbine
192    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  j_smear   !< index defining the area for the smearing of the forces (y-direction)
193    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  k_hub     !< index belonging to hub height
194    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  k_smear   !< index defining the area for the smearing of the forces (z-direction)
195    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nrings    !< number of rings per turbine
196    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nsegs_total !< total number of segments per turbine
[1819]197
[1912]198    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nsegs   !< number of segments per ring and turbine
[1819]199
[1912]200!
201!-  parameters for the smearing from the rotor to the cartesian grid   
[1864]202    REAL(wp) ::  pol_a            !< parameter for the polynomial smearing fct
203    REAL(wp) ::  pol_b            !< parameter for the polynomial smearing fct
[1912]204    REAL(wp) ::  delta_t_factor   !<
[1864]205    REAL(wp) ::  eps_factor       !< 
[1839]206    REAL(wp) ::  eps_min          !<
207    REAL(wp) ::  eps_min2         !<
208    REAL(wp) ::  sqrt_arg         !<
[1819]209
[1839]210!
211!-- Variables for the calculation of lift and drag coefficients
[1912]212    REAL(wp), DIMENSION(:), ALLOCATABLE  ::  ard     !<
213    REAL(wp), DIMENSION(:), ALLOCATABLE  ::  crd     !<
214    REAL(wp), DIMENSION(:), ALLOCATABLE  ::  delta_r !< radial segment length
215    REAL(wp), DIMENSION(:), ALLOCATABLE  ::  lrd     !<
[1864]216   
[1912]217    REAL(wp) ::  accu_cl_cd_tab = 0.1_wp  !< Accuracy of the interpolation of
218                                          !< the lift and drag coeff [deg]
[1819]219
[1912]220    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: turb_cd_tab   !< table of the blade drag coefficient
221    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: turb_cl_tab   !< table of the blade lift coefficient
[1819]222
[1912]223    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  nac_cd_surf  !< 3d field of the tower drag coefficient
224    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  tow_cd_surf  !< 3d field of the nacelle drag coefficient
[1819]225
226!
227!-- Variables for the calculation of the forces
[1912]228     
[1839]229    REAL(wp) ::  cur_r                       !<
[1912]230    REAL(wp) ::  delta_t                     !<  tangential segment length
[1839]231    REAL(wp) ::  phi_rotor                   !<
232    REAL(wp) ::  pre_factor                  !< 
233    REAL(wp) ::  torque_seg                  !<
234    REAL(wp) ::  u_int_l                     !<
235    REAL(wp) ::  u_int_u                     !<
236    REAL(wp) ::  u_rot                       !<
237    REAL(wp) ::  v_int_l                     !<
238    REAL(wp) ::  v_int_u                     !<
239    REAL(wp) ::  w_int_l                     !<
240    REAL(wp) ::  w_int_u                     !<
[1912]241!
242!-  Tendencies from the nacelle and tower thrust
243    REAL(wp) ::  tend_nac_x = 0.0_wp  !<
244    REAL(wp) ::  tend_tow_x = 0.0_wp  !<
245    REAL(wp) ::  tend_nac_y = 0.0_wp  !<
246    REAL(wp) ::  tend_tow_y = 0.0_wp  !<
[1819]247
[1912]248    REAL(wp), DIMENSION(:), ALLOCATABLE ::  alpha_attack !<
249    REAL(wp), DIMENSION(:), ALLOCATABLE ::  chord        !<
250    REAL(wp), DIMENSION(:), ALLOCATABLE ::  omega_gen    !< curr. generator speed
251    REAL(wp), DIMENSION(:), ALLOCATABLE ::  phi_rel      !<
252    REAL(wp), DIMENSION(:), ALLOCATABLE ::  pitch_add_old!<
253    REAL(wp), DIMENSION(:), ALLOCATABLE ::  torque_total !<
254    REAL(wp), DIMENSION(:), ALLOCATABLE ::  thrust_rotor !<
255    REAL(wp), DIMENSION(:), ALLOCATABLE ::  turb_cl      !<
256    REAL(wp), DIMENSION(:), ALLOCATABLE ::  turb_cd      !<
257    REAL(wp), DIMENSION(:), ALLOCATABLE ::  vrel         !<
258    REAL(wp), DIMENSION(:), ALLOCATABLE ::  vtheta       !<
259
260    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rbx, rby, rbz     !< coordinates of the blade elements
261    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rotx, roty, rotz  !< normal vectors to the rotor coordinates
262
263!
264!-  Fields for the interpolation of velocities on the rotor grid
265    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_int       !<
266    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_int_1_l   !<
267    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_int       !<
268    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_int_1_l   !<
269    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_int       !<
270    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_int_1_l   !<
271   
272!
273!-  rotor tendencies on the segments
274    REAL(wp), DIMENSION(:), ALLOCATABLE :: thrust_seg   !<
275    REAL(wp), DIMENSION(:), ALLOCATABLE :: torque_seg_y !<
276    REAL(wp), DIMENSION(:), ALLOCATABLE :: torque_seg_z !<   
277
278!
279!-  rotor tendencies on the rings
280    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  thrust_ring       !<
281    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  torque_ring_y     !<
282    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  torque_ring_z     !<
283   
284!
285!-  rotor tendencies on rotor grids for all turbines
286    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  thrust      !<
287    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  torque_y    !<
288    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  torque_z    !<
289
290!
291!-  rotor tendencies on coordinate grid
292    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rot_tend_x  !<
293    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rot_tend_y  !<
294    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rot_tend_z  !<
295!   
296!-  variables for the rotation of the rotor coordinates       
297    REAL(wp), DIMENSION(1:100,1:3,1:3) ::  rot_coord_trans  !< matrix for rotation of rotor coordinates
298   
[1839]299    REAL(wp), DIMENSION(1:3) ::  rot_eigen_rad   !<
300    REAL(wp), DIMENSION(1:3) ::  rot_eigen_azi   !<
301    REAL(wp), DIMENSION(1:3) ::  rot_eigen_nor   !<
302    REAL(wp), DIMENSION(1:3) ::  re              !<
303    REAL(wp), DIMENSION(1:3) ::  rea             !<
304    REAL(wp), DIMENSION(1:3) ::  ren             !<
305    REAL(wp), DIMENSION(1:3) ::  rote            !<
306    REAL(wp), DIMENSION(1:3) ::  rota            !<
307    REAL(wp), DIMENSION(1:3) ::  rotn            !<
[1819]308
[1839]309!
310!-- Fixed variables for the speed controller
[1819]311
[1912]312    LOGICAL  ::  start_up = .TRUE.   !<
[1864]313   
[1912]314    REAL(wp) ::  Fcorner             !< corner freq for the controller low pass filter
315    REAL(wp) ::  min_reg25           !< min region 2.5
316    REAL(wp) ::  om_rate             !< rotor speed change
317    REAL(wp) ::  slope15             !< slope in region 1.5
318    REAL(wp) ::  slope25             !< slope in region 2.5
319    REAL(wp) ::  trq_rate            !< torque change
320    REAL(wp) ::  vs_sysp             !<
321    REAL(wp) ::  lp_coeff            !< coeff for the controller low pass filter
[1864]322   
323    REAL(wp), DIMENSION(:), ALLOCATABLE :: omega_gen_old   !< last gen. speed
324    REAL(wp), DIMENSION(:), ALLOCATABLE :: omega_gen_f     !< filtered gen. sp
325    REAL(wp), DIMENSION(:), ALLOCATABLE :: omega_gen_f_old !< last filt. gen. sp
326    REAL(wp), DIMENSION(:), ALLOCATABLE :: torque_gen      !< generator torque
327    REAL(wp), DIMENSION(:), ALLOCATABLE :: torque_gen_old  !< last gen. torque
[1819]328
[1864]329    REAL(wp), DIMENSION(100) :: omega_rot_l = 0.0_wp !< local rot speed [rad/s]
[1839]330!
331!-- Fixed variables for the yaw controller
[1819]332
[1912]333    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  yawdir           !< direction to yaw
334    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  phi_yaw_l        !< local (cpu) yaw angle
335    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  wd30_l           !< local (cpu) long running avg of the wd
336    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  wd2_l            !< local (cpu) short running avg of the wd
337    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  wdir             !< wind direction at hub
338    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  u_inflow         !< wind speed at hub
339    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  wdir_l           !<
340    REAL(wp), DIMENSION(:)  , ALLOCATABLE ::  u_inflow_l       !<
341    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wd30             !<
342    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wd2              !<
343    LOGICAL,  DIMENSION(1:100)            ::  doyaw = .FALSE.  !<
344    INTEGER(iwp)                          ::  WDLON            !<
345    INTEGER(iwp)                          ::  WDSHO            !<
[1819]346
347
[1839]348    SAVE
[1819]349
[1839]350
351    INTERFACE wtm_parin
352       MODULE PROCEDURE wtm_parin
353    END INTERFACE wtm_parin
[1912]354   
355    INTERFACE wtm_check_parameters
356       MODULE PROCEDURE wtm_check_parameters
357    END INTERFACE wtm_check_parameters
[1839]358       
359    INTERFACE wtm_init_arrays
360       MODULE PROCEDURE wtm_init_arrays
361    END INTERFACE wtm_init_arrays
362
363    INTERFACE wtm_init
364       MODULE PROCEDURE wtm_init
365    END INTERFACE wtm_init
366   
[1912]367    INTERFACE wtm_read_blade_tables
368       MODULE PROCEDURE wtm_read_blade_tables
369    END INTERFACE wtm_read_blade_tables
[1839]370           
[1819]371    INTERFACE wtm_forces
372       MODULE PROCEDURE wtm_forces
[1839]373       MODULE PROCEDURE wtm_yawcontrol
[1819]374    END INTERFACE wtm_forces
[1864]375   
376    INTERFACE wtm_rotate_rotor
377       MODULE PROCEDURE wtm_rotate_rotor
378    END INTERFACE wtm_rotate_rotor
379   
380    INTERFACE wtm_speed_control
381       MODULE PROCEDURE wtm_init_speed_control
382       MODULE PROCEDURE wtm_speed_control
383    END INTERFACE wtm_speed_control
[1819]384
385    INTERFACE wtm_tendencies
386       MODULE PROCEDURE wtm_tendencies
387       MODULE PROCEDURE wtm_tendencies_ij
388    END INTERFACE wtm_tendencies
[1864]389   
390   
[1912]391    PUBLIC wtm_check_parameters, wtm_forces, wtm_init, wtm_init_arrays,        &
392           wtm_parin, wtm_tendencies, wtm_tendencies_ij, wind_turbine 
[1819]393
394 CONTAINS
395
396
397!------------------------------------------------------------------------------!
398! Description:
399! ------------
[1839]400!> Parin for &wind_turbine_par for wind turbine model
[1819]401!------------------------------------------------------------------------------!
[1839]402    SUBROUTINE wtm_parin
[1819]403
404
405       IMPLICIT NONE
[1912]406       
407       INTEGER(iwp) ::  ierrn       !<
[1819]408
[1839]409       CHARACTER (LEN=80) ::  line  !< dummy string that contains the current line of the parameter file
410
411       NAMELIST /wind_turbine_par/   air_dens, dtow, gear_eff, gear_ratio,     &
[1864]412                                  gen_eff, inertia_gen, inertia_rot, max_miss, &
[1839]413                                  max_torque_gen, max_trq_rate, min_miss,      &
[1864]414                                  min_reg15, min_reg2, nairfoils, nturbines,   &
[1912]415                                  omega_rot, phi_yaw, pitch_add, pitch_control,&
[1839]416                                  rated_genspeed, rated_power, rcx, rcy, rcz,  &
417                                  rnac, rr, segment_length, segment_width,     &
[1864]418                                  slope2, speed_control, tilt, time_turbine_on,&
[1912]419                                  turb_cd_nacelle, turb_cd_tower,              &
[2152]420                                  yaw_control, yaw_speed, tl_cor
[1839]421
[1819]422!
[1839]423!--    Try to find wind turbine model package
424       REWIND ( 11 )
425       line = ' '
426       DO  WHILE ( INDEX( line, '&wind_turbine_par' ) == 0 )
427          READ ( 11, '(A)', END=10 )  line
428       ENDDO
429       BACKSPACE ( 11 )
430
431!
432!--    Read user-defined namelist
[1912]433       READ ( 11, wind_turbine_par, IOSTAT=ierrn )
434       
435       IF ( ierrn < 0 )  THEN
436          message_string = 'errors in \$wind_turbine_par'
[2322]437          CALL message( 'wtm_parin', 'PA0460', 1, 2, 0, 6, 0 )
[1912]438       ENDIF
439       
[1839]440!
441!--    Set flag that indicates that the wind turbine model is switched on
442       wind_turbine = .TRUE.
443
[1912]444 10    CONTINUE   ! TBD Change from continue, mit ierrn machen
[1839]445
446
447    END SUBROUTINE wtm_parin
448
[1912]449    SUBROUTINE wtm_check_parameters
450
451   
452       IMPLICIT NONE
453   
454       IF ( ( .NOT.speed_control ) .AND. pitch_control )  THEN
455          message_string = 'pitch_control = .TRUE. requires '//                &
456                           'speed_control = .TRUE.'
[2322]457          CALL message( 'wtm_check_parameters', 'PA0461', 1, 2, 0, 6, 0 )
[1912]458       ENDIF
459       
[2322]460       IF ( ANY( omega_rot(1:nturbines) < 0.0 ) )  THEN
461          message_string = 'omega_rot < 0.0, Please set omega_rot to '     // &
[1912]462                           'a value larger than zero'
[2322]463          CALL message( 'wtm_check_parameters', 'PA0462', 1, 2, 0, 6, 0 )
[1912]464       ENDIF
465       
466       
467       IF ( ANY( rcx(1:nturbines) == 9999999.9_wp ) .OR.                       &
468            ANY( rcy(1:nturbines) == 9999999.9_wp ) .OR.                       &
469            ANY( rcz(1:nturbines) == 9999999.9_wp ) )  THEN
470         
471          message_string = 'rcx, rcy, rcz '                                 // &
472                           'have to be given for each turbine.'         
[2322]473          CALL message( 'wtm_check_parameters', 'PA0463', 1, 2, 0, 6, 0 )         
[1912]474         
475       ENDIF
476
477 
478    END SUBROUTINE wtm_check_parameters 
479   
480                                       
[1839]481!------------------------------------------------------------------------------!
482! Description:
483! ------------
484!> Allocate wind turbine model arrays
485!------------------------------------------------------------------------------!
486    SUBROUTINE wtm_init_arrays
487
488
489       IMPLICIT NONE
490
[1864]491       REAL(wp) ::  delta_r_factor   !<
492       REAL(wp) ::  delta_r_init     !<
493
[1839]494!
495!--    To be able to allocate arrays with dimension of rotor rings and segments,
[1819]496!--    the maximum possible numbers of rings and segments have to be calculated:
497
498       ALLOCATE( nrings(1:nturbines) )
499       ALLOCATE( delta_r(1:nturbines) )
500
501       nrings(:)  = 0
502       delta_r(:) = 0.0_wp
503
504!
505!--    Thickness (radial) of each ring and length (tangential) of each segment:
506       delta_r_factor = segment_width
507       delta_t_factor = segment_length
508       delta_r_init   = delta_r_factor * MIN( dx, dy, dz)
509       delta_t        = delta_t_factor * MIN( dx, dy, dz)
510
511       DO inot = 1, nturbines
512!
513!--       Determine number of rings:
514          nrings(inot) = NINT( rr(inot) / delta_r_init )
515
516          delta_r(inot) = rr(inot) / nrings(inot)
517
518       ENDDO
519
520       nrings_max = MAXVAL(nrings)
521
522       ALLOCATE( nsegs(1:nrings_max,1:nturbines) )
523       ALLOCATE( nsegs_total(1:nturbines) )
524
525       nsegs(:,:)     = 0
526       nsegs_total(:) = 0
527
528
529       DO inot = 1, nturbines
530          DO ring = 1, nrings(inot)
531!
532!--          Determine number of segments for each ring:
[1839]533             nsegs(ring,inot) = MAX( 8, CEILING( delta_r_factor * pi *         &
534                                                 ( 2.0_wp * ring - 1.0_wp ) /  &
535                                                 delta_t_factor ) )
[1819]536          ENDDO
537!
538!--       Total sum of all rotor segments:
[1839]539          nsegs_total(inot) = SUM( nsegs(:,inot) )
[1819]540
541       ENDDO
542
543!
544!--    Maximum number of segments per ring:
545       nsegs_max = MAXVAL(nsegs)
546
[1864]547!!
548!!--    TODO: Folgendes im Header ausgeben!
549!       IF ( myid == 0 )  THEN
550!          PRINT*, 'nrings(1) = ', nrings(1)
551!          PRINT*, '--------------------------------------------------'
552!          PRINT*, 'nsegs(:,1) = ', nsegs(:,1)
553!          PRINT*, '--------------------------------------------------'
554!          PRINT*, 'nrings_max = ', nrings_max
555!          PRINT*, 'nsegs_max = ', nsegs_max
556!          PRINT*, 'nsegs_total(1) = ', nsegs_total(1)
557!       ENDIF
[1819]558
559
560!
561!--    Allocate 1D arrays (dimension = number of turbines)
562       ALLOCATE( i_hub(1:nturbines) )
563       ALLOCATE( i_smear(1:nturbines) )
564       ALLOCATE( j_hub(1:nturbines) )
565       ALLOCATE( j_smear(1:nturbines) )
566       ALLOCATE( k_hub(1:nturbines) )
567       ALLOCATE( k_smear(1:nturbines) )
568       ALLOCATE( torque_total(1:nturbines) )
[1912]569       ALLOCATE( thrust_rotor(1:nturbines) )
[1819]570
571!
572!--    Allocation of the 1D arrays for speed pitch_control
573       ALLOCATE( omega_gen(1:nturbines) )
[1864]574       ALLOCATE( omega_gen_old(1:nturbines) )
[1819]575       ALLOCATE( omega_gen_f(1:nturbines) )
576       ALLOCATE( omega_gen_f_old(1:nturbines) )
[1912]577       ALLOCATE( pitch_add_old(1:nturbines) )
[1819]578       ALLOCATE( torque_gen(1:nturbines) )
579       ALLOCATE( torque_gen_old(1:nturbines) )
580
581!
582!--    Allocation of the 1D arrays for yaw control
583       ALLOCATE( yawdir(1:nturbines) )
584       ALLOCATE( u_inflow(1:nturbines) )
585       ALLOCATE( wdir(1:nturbines) )
586       ALLOCATE( u_inflow_l(1:nturbines) )
587       ALLOCATE( wdir_l(1:nturbines) )
588       ALLOCATE( phi_yaw_l(1:nturbines) )
589       
590!
591!--    Allocate 1D arrays (dimension = number of rotor segments)
592       ALLOCATE( alpha_attack(1:nsegs_max) )
593       ALLOCATE( chord(1:nsegs_max) )
594       ALLOCATE( phi_rel(1:nsegs_max) )
595       ALLOCATE( thrust_seg(1:nsegs_max) )
596       ALLOCATE( torque_seg_y(1:nsegs_max) )
597       ALLOCATE( torque_seg_z(1:nsegs_max) )
[1839]598       ALLOCATE( turb_cd(1:nsegs_max) )
599       ALLOCATE( turb_cl(1:nsegs_max) )
[1819]600       ALLOCATE( vrel(1:nsegs_max) )
601       ALLOCATE( vtheta(1:nsegs_max) )
602
603!
604!--    Allocate 2D arrays (dimension = number of rotor rings and segments)
605       ALLOCATE( rbx(1:nrings_max,1:nsegs_max) )
606       ALLOCATE( rby(1:nrings_max,1:nsegs_max) )
607       ALLOCATE( rbz(1:nrings_max,1:nsegs_max) )
608       ALLOCATE( thrust_ring(1:nrings_max,1:nsegs_max) )
609       ALLOCATE( torque_ring_y(1:nrings_max,1:nsegs_max) )
610       ALLOCATE( torque_ring_z(1:nrings_max,1:nsegs_max) )
611
612!
613!--    Allocate additional 2D arrays
614       ALLOCATE( rotx(1:nturbines,1:3) )
615       ALLOCATE( roty(1:nturbines,1:3) )
616       ALLOCATE( rotz(1:nturbines,1:3) )
617
618!
619!--    Allocate 3D arrays (dimension = number of grid points)
[1912]620       ALLOCATE( nac_cd_surf(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
621       ALLOCATE( rot_tend_x(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
622       ALLOCATE( rot_tend_y(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
623       ALLOCATE( rot_tend_z(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1819]624       ALLOCATE( thrust(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
625       ALLOCATE( torque_y(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
626       ALLOCATE( torque_z(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1912]627       ALLOCATE( tow_cd_surf(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1819]628
629!
630!--    Allocate additional 3D arrays
631       ALLOCATE( u_int(1:nturbines,1:nrings_max,1:nsegs_max) )
632       ALLOCATE( u_int_1_l(1:nturbines,1:nrings_max,1:nsegs_max) )
633       ALLOCATE( v_int(1:nturbines,1:nrings_max,1:nsegs_max) )
634       ALLOCATE( v_int_1_l(1:nturbines,1:nrings_max,1:nsegs_max) )
635       ALLOCATE( w_int(1:nturbines,1:nrings_max,1:nsegs_max) )
636       ALLOCATE( w_int_1_l(1:nturbines,1:nrings_max,1:nsegs_max) )
637
638!
639!--    All of the arrays are initialized with a value of zero:
640       i_hub(:)                 = 0
641       i_smear(:)               = 0
642       j_hub(:)                 = 0
643       j_smear(:)               = 0
644       k_hub(:)                 = 0
645       k_smear(:)               = 0
[1912]646       
[1819]647       torque_total(:)          = 0.0_wp
[1912]648       thrust_rotor(:)          = 0.0_wp
[1819]649
650       omega_gen(:)             = 0.0_wp
[1864]651       omega_gen_old(:)         = 0.0_wp
[1819]652       omega_gen_f(:)           = 0.0_wp
653       omega_gen_f_old(:)       = 0.0_wp
[1912]654       pitch_add_old(:)         = 0.0_wp
[1819]655       torque_gen(:)            = 0.0_wp
656       torque_gen_old(:)        = 0.0_wp
657       
658       yawdir(:)                = 0.0_wp
659       wdir(:)                  = 0.0_wp
660       u_inflow(:)              = 0.0_wp
661
662!
663!--    Allocate 1D arrays (dimension = number of rotor segments)
664       alpha_attack(:)          = 0.0_wp
665       chord(:)                 = 0.0_wp
666       phi_rel(:)               = 0.0_wp
667       thrust_seg(:)            = 0.0_wp
668       torque_seg_y(:)          = 0.0_wp
669       torque_seg_z(:)          = 0.0_wp
[1864]670       turb_cd(:)               = 0.0_wp
671       turb_cl(:)               = 0.0_wp
[1819]672       vrel(:)                  = 0.0_wp
673       vtheta(:)                = 0.0_wp
674
675       rbx(:,:)                 = 0.0_wp
676       rby(:,:)                 = 0.0_wp
677       rbz(:,:)                 = 0.0_wp
678       thrust_ring(:,:)         = 0.0_wp
679       torque_ring_y(:,:)       = 0.0_wp
680       torque_ring_z(:,:)       = 0.0_wp
681
682       rotx(:,:)                = 0.0_wp
683       roty(:,:)                = 0.0_wp
684       rotz(:,:)                = 0.0_wp
685
[1912]686       nac_cd_surf(:,:,:)       = 0.0_wp
687       rot_tend_x(:,:,:)        = 0.0_wp
688       rot_tend_y(:,:,:)        = 0.0_wp
689       rot_tend_z(:,:,:)        = 0.0_wp
[1819]690       thrust(:,:,:)            = 0.0_wp
691       torque_y(:,:,:)          = 0.0_wp
692       torque_z(:,:,:)          = 0.0_wp
[1912]693       tow_cd_surf(:,:,:)       = 0.0_wp
[1819]694
695       u_int(:,:,:)             = 0.0_wp
696       u_int_1_l(:,:,:)         = 0.0_wp
697       v_int(:,:,:)             = 0.0_wp
698       v_int_1_l(:,:,:)         = 0.0_wp
699       w_int(:,:,:)             = 0.0_wp
700       w_int_1_l(:,:,:)         = 0.0_wp
701
702
[1839]703    END SUBROUTINE wtm_init_arrays
[1819]704
705
706!------------------------------------------------------------------------------!
707! Description:
708! ------------
[1839]709!> Initialization of the wind turbine model
[1819]710!------------------------------------------------------------------------------!
[1839]711    SUBROUTINE wtm_init
[1819]712
[1839]713   
[1819]714       IMPLICIT NONE
715
716       INTEGER(iwp) ::  i  !< running index
717       INTEGER(iwp) ::  j  !< running index
718       INTEGER(iwp) ::  k  !< running index
[1864]719       
[1819]720!
[1864]721!--    Help variables for the smearing function       
722       REAL(wp) ::  eps_kernel       !<       
723       
[1839]724!
[1864]725!--    Help variables for calculation of the tower drag       
726       INTEGER(iwp) ::  tower_n      !<
727       INTEGER(iwp) ::  tower_s      !<
[1839]728!
[1864]729!--    Help variables for the calulaction of the nacelle drag
730       INTEGER(iwp) ::  i_ip         !<
731       INTEGER(iwp) ::  i_ipg        !<
732       
733       REAL(wp) ::  yvalue               
734       REAL(wp) ::  dy_int           !<
735       REAL(wp) ::  dz_int           !<
736       
737       REAL(wp), DIMENSION(:,:), ALLOCATABLE :: circle_points  !<
738             
739       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: index_nacb       !<
740       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: index_nacl       !<
741       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: index_nacr       !<
742       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: index_nact       !<
743       
744       ALLOCATE( index_nacb(1:nturbines) )
745       ALLOCATE( index_nacl(1:nturbines) )
746       ALLOCATE( index_nacr(1:nturbines) )
747       ALLOCATE( index_nact(1:nturbines) )
[1819]748
749
[1864]750       IF ( speed_control)  THEN
751       
752          CALL wtm_speed_control
753
754       ENDIF
755
[1839]756!
[1819]757!------------------------------------------------------------------------------!
[1839]758!--    Calculation of parameters for the regularization kernel
759!--    (smearing of the forces)
[1819]760!------------------------------------------------------------------------------!
761!
[1839]762!--    In the following, some of the required parameters for the smearing will
763!--    be calculated:
[1819]764
[1839]765!--    The kernel is set equal to twice the grid spacing which has turned out to
766!--    be a reasonable value (see e.g. Troldborg et al. (2013), Wind Energy,
[1819]767!--    DOI: 10.1002/we.1608):
768       eps_kernel = 2.0_wp * dx
769!
[1839]770!--    The zero point (eps_min) of the polynomial function must be the following
771!--    if the integral of the polynomial function (for values < eps_min) shall
772!--    be equal to the integral of the Gaussian function used before:
773       eps_min = ( 105.0_wp / 32.0_wp )**( 1.0_wp / 3.0_wp ) *                 &
774                 pi**( 1.0_wp / 6.0_wp ) * eps_kernel
[1819]775!
776!--    Square of eps_min:
777       eps_min2 = eps_min**2
778!
779!--    Parameters in the polynomial function:
[1864]780       pol_a = 1.0_wp / eps_min**4
781       pol_b = 2.0_wp / eps_min**2
[1819]782!
[1839]783!--    Normalization factor which is the inverse of the integral of the smearing
784!--    function:
785       eps_factor = 105.0_wp / ( 32.0_wp * pi * eps_min**3 )
786       
[1864]787!--    Change tilt angle to rad:
788       tilt = tilt * pi / 180.0_wp
789     
[1819]790!
[1864]791!--    Change yaw angle to rad:
792       phi_yaw(:) = phi_yaw(:) * pi / 180.0_wp
[1819]793
[1864]794
[1819]795       DO inot = 1, nturbines
796!
[1864]797!--       Rotate the rotor coordinates in case yaw and tilt are defined
798          CALL wtm_rotate_rotor( inot )
[1819]799         
800!
801!--       Determine the indices of the hub height
802          i_hub(inot) = INT(   rcx(inot)                 / dx )
803          j_hub(inot) = INT( ( rcy(inot) + 0.5_wp * dy ) / dy )
804          k_hub(inot) = INT( ( rcz(inot) + 0.5_wp * dz ) / dz )
805
806!
807!--       Determining the area to which the smearing of the forces is applied.
[1839]808!--       As smearing now is effectively applied only for distances smaller than
809!--       eps_min, the smearing area can be further limited and regarded as a
810!--       function of eps_min:
[1819]811          i_smear(inot) = CEILING( ( rr(inot) + eps_min ) / dx )
812          j_smear(inot) = CEILING( ( rr(inot) + eps_min ) / dy )
813          k_smear(inot) = CEILING( ( rr(inot) + eps_min ) / dz )
[1864]814       
[1819]815       ENDDO
816
817!
818!------------------------------------------------------------------------------!
[1839]819!--    Determine the area within each grid cell that overlaps with the area
820!--    of the nacelle and the tower (needed for calculation of the forces)
[1819]821!------------------------------------------------------------------------------!
822!
823!--    Note: so far this is only a 2D version, in that the mean flow is
824!--    perpendicular to the rotor area.
825
826!
827!--    Allocation of the array containing information on the intersection points
828!--    between rotor disk and the numerical grid:
829       upper_end = ( ny + 1 ) * 10000 
830
831       ALLOCATE( circle_points(1:2,1:upper_end) )
[1839]832       
833       circle_points(:,:) = 0.0_wp
[1819]834
[1839]835       
836       DO inot = 1, nturbines                     ! loop over number of turbines
[1819]837!
[1839]838!--       Determine the grid index (u-grid) that corresponds to the location of
839!--       the rotor center (reduces the amount of calculations in the case that
840!--       the mean flow is perpendicular to the rotor area):
[1819]841          i = i_hub(inot)
842
843!
[1839]844!--       Determine the left and the right edge of the nacelle (corresponding
845!--       grid point indices):
[1819]846          index_nacl(inot) = INT( ( rcy(inot) - rnac(inot) + 0.5_wp * dy ) / dy )
847          index_nacr(inot) = INT( ( rcy(inot) + rnac(inot) + 0.5_wp * dy ) / dy )
848!
[1839]849!--       Determine the bottom and the top edge of the nacelle (corresponding
850!--       grid point indices).The grid point index has to be increased by 1, as
851!--       the first level for the u-component (index 0) is situated below the
852!--       surface. All points between z=0 and z=dz/s would already be contained
853!--       in grid box 1.
[1819]854          index_nacb(inot) = INT( ( rcz(inot) - rnac(inot) ) / dz ) + 1
855          index_nact(inot) = INT( ( rcz(inot) + rnac(inot) ) / dz ) + 1
856
857!
858!--       Determine the indices of the grid boxes containing the left and
[1864]859!--       the right boundaries of the tower:
[1819]860          tower_n = ( rcy(inot) + 0.5_wp * dtow(inot) - 0.5_wp * dy ) / dy
861          tower_s = ( rcy(inot) - 0.5_wp * dtow(inot) - 0.5_wp * dy ) / dy
862
863!
864!--       Determine the fraction of the grid box area overlapping with the tower
[1864]865!--       area and multiply it with the drag of the tower:
[1839]866          IF ( ( nxlg <= i )  .AND.  ( nxrg >= i ) )  THEN
[1819]867
868             DO  j = nys, nyn
869!
[1839]870!--             Loop from south to north boundary of tower
871                IF ( ( j >= tower_s )  .AND.  ( j <= tower_n ) )  THEN
872
[1819]873                   DO  k = nzb, nzt
874
875                      IF ( k == k_hub(inot) )  THEN
876                         IF ( tower_n - tower_s >= 1 )  THEN
[1839]877!
[1819]878!--                      leftmost and rightmost grid box:
879                            IF ( j == tower_s )  THEN
[1912]880                               tow_cd_surf(k,j,i) = ( rcz(inot) -              &
881                                    ( k_hub(inot) * dz - 0.5_wp * dz ) )  *    & ! extension in z-direction
882                                  ( ( tower_s + 1.0_wp + 0.5_wp ) * dy    -    &
883                                    ( rcy(inot) - 0.5_wp * dtow(inot) ) ) *    & ! extension in y-direction
884                                  turb_cd_tower(inot)
[1819]885                            ELSEIF ( j == tower_n )  THEN
[1912]886                               tow_cd_surf(k,j,i) = ( rcz(inot)            -   &
887                                    ( k_hub(inot) * dz - 0.5_wp * dz ) )  *    & ! extension in z-direction
888                                  ( ( rcy(inot) + 0.5_wp * dtow(inot) )   -    &
889                                    ( tower_n + 0.5_wp ) * dy )           *    & ! extension in y-direction
890                                  turb_cd_tower(inot)
[1819]891!
892!--                         grid boxes inbetween
[1912]893!--                         (where tow_cd_surf = grid box area):
[1819]894                            ELSE
[1912]895                               tow_cd_surf(k,j,i) = ( rcz(inot) -              &
896                                    ( k_hub(inot) * dz - 0.5_wp * dz ) )  *    &
897                                    dy * turb_cd_tower(inot)
[1819]898                            ENDIF
899!
900!--                      tower lies completely within one grid box:
901                         ELSE
[1912]902                            tow_cd_surf(k,j,i) = ( rcz(inot)                 - &
[1839]903                                       ( k_hub(inot) * dz - 0.5_wp * dz ) ) *  &
[1912]904                                       dtow(inot) * turb_cd_tower(inot)
[1819]905                         ENDIF
[1839]906!
907!--                  In case that k is smaller than k_hub the following actions
908!--                  are carried out:
[1819]909                      ELSEIF ( k < k_hub(inot) )  THEN
910                     
[1839]911                         IF ( ( tower_n - tower_s ) >= 1 )  THEN
912!
[1819]913!--                         leftmost and rightmost grid box:
914                            IF ( j == tower_s )  THEN                         
[1912]915                               tow_cd_surf(k,j,i) = dz * (                     &
[1839]916                                      ( tower_s + 1 + 0.5_wp ) * dy         -  &
917                                      ( rcy(inot) - 0.5_wp * dtow(inot) )      &
[1912]918                                                        ) * turb_cd_tower(inot)
[1819]919                            ELSEIF ( j == tower_n )  THEN
[1912]920                               tow_cd_surf(k,j,i) = dz * (                     &
[1839]921                                      ( rcy(inot) + 0.5_wp * dtow(inot) )   -  &
922                                      ( tower_n + 0.5_wp ) * dy                &
[1912]923                                                         ) * turb_cd_tower(inot)
[1839]924!
925!--                         grid boxes inbetween
[1912]926!--                         (where tow_cd_surf = grid box area):
[1819]927                            ELSE
[1912]928                               tow_cd_surf(k,j,i) = dz * dy * turb_cd_tower(inot)
[1819]929                            ENDIF
[1839]930!
[1819]931!--                         tower lies completely within one grid box:
932                         ELSE
[1912]933                            tow_cd_surf(k,j,i) = dz * dtow(inot) *             &
934                                                turb_cd_tower(inot)
[1839]935                         ENDIF ! end if larger than grid box
936
937                      ENDIF    ! end if k == k_hub
938
939                   ENDDO       ! end loop over k
940
941                ENDIF          ! end if inside north and south boundary of tower
942
943             ENDDO             ! end loop over j
944
945          ENDIF                ! end if hub inside domain + ghostpoints
[1819]946       
[1839]947         
[1912]948          CALL exchange_horiz( tow_cd_surf, nbgp )
[1819]949
[1839]950!
[1864]951!--       Calculation of the nacelle area
952!--       CAUTION: Currently disabled due to segmentation faults on the FLOW HPC
953!--                cluster (Oldenburg)
954!!
955!!--       Tabulate the points on the circle that are required in the following for
956!!--       the calculation of the Riemann integral (node points; they are called
957!!--       circle_points in the following):
958!
959!          dy_int = dy / 10000.0_wp
960!
961!          DO  i_ip = 1, upper_end
962!             yvalue   = dy_int * ( i_ip - 0.5_wp ) + 0.5_wp * dy           !<--- segmentation fault
963!             sqrt_arg = rnac(inot)**2 - ( yvalue - rcy(inot) )**2          !<--- segmentation fault
964!             IF ( sqrt_arg >= 0.0_wp )  THEN
965!!
966!!--             bottom intersection point
967!                circle_points(1,i_ip) = rcz(inot) - SQRT( sqrt_arg )
968!!
969!!--             top intersection point
970!                circle_points(2,i_ip) = rcz(inot) + SQRT( sqrt_arg )       !<--- segmentation fault
971!             ELSE
972!                circle_points(:,i_ip) = -111111                            !<--- segmentation fault
973!             ENDIF
974!          ENDDO
975!
976!
977!          DO  j = nys, nyn
978!!
979!!--          In case that the grid box is located completely outside the nacelle
980!!--          (y) it can automatically be stated that there is no overlap between
981!!--          the grid box and the nacelle and consequently we can set
[1912]982!!--          nac_cd_surf(:,j,i) = 0.0:
[1864]983!             IF ( ( j >= index_nacl(inot) )  .AND.  ( j <= index_nacr(inot) ) )  THEN
984!                DO  k = nzb+1, nzt
985!!
986!!--                In case that the grid box is located completely outside the
987!!--                nacelle (z) it can automatically be stated that there is no
988!!--                overlap between the grid box and the nacelle and consequently
[1912]989!!--                we can set nac_cd_surf(k,j,i) = 0.0:
[1864]990!                   IF ( ( k >= index_nacb(inot) )  .OR.                           &
991!                        ( k <= index_nact(inot) ) )  THEN
992!!
993!!--                   For all other cases Riemann integrals are calculated.
994!!--                   Here, the points on the circle that have been determined
995!!--                   above are used in order to calculate the overlap between the
996!!--                   gridbox and the nacelle area (area approached by 10000
997!!--                   rectangulars dz_int * dy_int):
998!                      DO  i_ipg = 1, 10000
999!                         dz_int = dz
1000!                         i_ip = j * 10000 + i_ipg
1001!!
1002!!--                      Determine the vertical extension dz_int of the circle
1003!!--                      within the current grid box:
1004!                         IF ( ( circle_points(2,i_ip) < zw(k) ) .AND.          &  !<--- segmentation fault
1005!                              ( circle_points(2,i_ip) >= zw(k-1) ) ) THEN
1006!                            dz_int = dz_int -                                  &  !<--- segmentation fault
1007!                                     ( zw(k) - circle_points(2,i_ip) )
1008!                         ENDIF
1009!                         IF ( ( circle_points(1,i_ip) <= zw(k) ) .AND.         &  !<--- segmentation fault
1010!                              ( circle_points(1,i_ip) > zw(k-1) ) ) THEN
1011!                            dz_int = dz_int -                                  &
1012!                                     ( circle_points(1,i_ip) - zw(k-1) )
1013!                         ENDIF
1014!                         IF ( zw(k-1) > circle_points(2,i_ip) ) THEN
1015!                            dz_int = 0.0_wp
1016!                         ENDIF
1017!                         IF ( zw(k) < circle_points(1,i_ip) ) THEN
1018!                            dz_int = 0.0_wp                     
1019!                         ENDIF
1020!                         IF ( ( nxlg <= i ) .AND. ( nxrg >= i ) ) THEN
[1912]1021!                            nac_cd_surf(k,j,i) = nac_cd_surf(k,j,i) +        &  !<--- segmentation fault
1022!                                                  dy_int * dz_int * turb_cd_nacelle(inot)
[1864]1023!                         ENDIF   
1024!                      ENDDO
1025!                   ENDIF
1026!                ENDDO
1027!             ENDIF
[1819]1028!
[1864]1029!          ENDDO
[1819]1030!       
[1912]1031!          CALL exchange_horiz( nac_cd_surf, nbgp )                                !<---  segmentation fault
[1819]1032
[1864]1033       ENDDO   ! end of loop over turbines
[1819]1034
[1912]1035       tow_cd_surf   = tow_cd_surf   / ( dx * dy * dz )      ! Normalize tower drag
1036       nac_cd_surf = nac_cd_surf / ( dx * dy * dz )      ! Normalize nacelle drag
[1819]1037
[1912]1038       CALL wtm_read_blade_tables
[1839]1039 
[1864]1040    END SUBROUTINE wtm_init
1041
1042
[1819]1043!------------------------------------------------------------------------------!
[1864]1044! Description:
1045! ------------
1046!> Read in layout of the rotor blade , the lift and drag tables
1047!> and the distribution of lift and drag tables along the blade
[1819]1048!------------------------------------------------------------------------------!
[1912]1049!
1050    SUBROUTINE wtm_read_blade_tables
[1819]1051
[1839]1052
[1864]1053       IMPLICIT NONE
1054
1055       INTEGER(iwp) ::  ii   !< running index
1056       INTEGER(iwp) ::  jj   !< running index
[1843]1057   
[1864]1058       INTEGER(iwp) ::  ierrn       !<
1059   
1060       CHARACTER(200) :: chmess     !< Read in string
[1839]1061
[1864]1062       INTEGER(iwp) ::  dlen        !< no. rows of local table
1063       INTEGER(iwp) ::  dlenbl      !< no. rows of cd, cl table
1064       INTEGER(iwp) ::  ialpha      !< table position of current alpha value
1065       INTEGER(iwp) ::  iialpha     !<
1066       INTEGER(iwp) ::  iir         !<
1067       INTEGER(iwp) ::  radres      !< radial resolution
1068       INTEGER(iwp) ::  t1          !< no. of airfoil
1069       INTEGER(iwp) ::  t2          !< no. of airfoil
1070       INTEGER(iwp) ::  trow        !<
1071       INTEGER(iwp) ::  dlenbl_int  !< no. rows of interpolated cd, cl tables
[1839]1072   
[1864]1073       REAL(wp) :: alpha_attack_i   !<
1074       REAL(wp) :: weight_a         !<
1075       REAL(wp) :: weight_b         !<
[1839]1076
[1864]1077       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: ttoint1    !<
1078       INTEGER(iwp), DIMENSION(:), ALLOCATABLE :: ttoint2    !<
[1839]1079   
[1864]1080       REAL(wp), DIMENSION(:), ALLOCATABLE :: turb_cd_sel1   !<
1081       REAL(wp), DIMENSION(:), ALLOCATABLE :: turb_cd_sel2   !<
1082       REAL(wp), DIMENSION(:), ALLOCATABLE :: turb_cl_sel1   !<
1083       REAL(wp), DIMENSION(:), ALLOCATABLE :: turb_cl_sel2   !<
1084       REAL(wp), DIMENSION(:), ALLOCATABLE :: read_cl_cd     !< read in var array
[1839]1085             
[1864]1086       REAL(wp), DIMENSION(:), ALLOCATABLE    :: alpha_attack_tab   !<
1087       REAL(wp), DIMENSION(:), ALLOCATABLE    :: trad1              !<
1088       REAL(wp), DIMENSION(:), ALLOCATABLE    :: trad2              !<         
1089       REAL(wp), DIMENSION(:,:), ALLOCATABLE  :: turb_cd_table      !<
1090       REAL(wp), DIMENSION(:,:), ALLOCATABLE  :: turb_cl_table      !<
[1839]1091                                         
[1864]1092       ALLOCATE ( read_cl_cd(1:2*nairfoils+1) )
[1839]1093
[1864]1094!
1095!--    Read in the distribution of lift and drag tables along the blade, the
1096!--    layout of the rotor blade and the lift and drag tables:
1097
1098       OPEN ( 201, FILE='WTM_DATA', STATUS='OLD', FORM='FORMATTED', IOSTAT=ierrn )
1099
1100       IF ( ierrn /= 0 )  THEN
1101          message_string = 'file WTM_DATA does not exist'
[2322]1102          CALL message( 'wtm_init', 'PA0464', 1, 2, 0, 6, 0 )
[1864]1103       ENDIF
1104!
1105!--    Read distribution table:
1106
1107       dlen = 0
1108
1109       READ ( 201, '(3/)' )
1110
1111       rloop3: DO
1112          READ ( 201, *, IOSTAT=ierrn ) chmess
1113          IF ( ierrn < 0  .OR.  chmess == '#'  .OR.  chmess == '')  EXIT rloop3
1114          dlen = dlen + 1
1115       ENDDO rloop3
1116
1117       ALLOCATE( trad1(1:dlen), trad2(1:dlen), ttoint1(1:dlen), ttoint2(1:dlen))
1118
1119       DO jj = 1,dlen+1
1120          BACKSPACE ( 201, IOSTAT=ierrn )
1121       ENDDO
1122
1123       DO jj = 1,dlen
1124          READ ( 201, * ) trad1(jj), trad2(jj), ttoint1(jj), ttoint2(jj)
1125       ENDDO
1126
1127!
1128!--    Read layout table:
1129
[1839]1130       dlen = 0 
[1864]1131
[1839]1132       READ ( 201, '(3/)')
[1864]1133
[1843]1134       rloop1: DO
[1864]1135          READ ( 201, *, IOSTAT=ierrn ) chmess
1136          IF ( ierrn < 0  .OR.  chmess == '#'  .OR.  chmess == '')  EXIT rloop1
1137          dlen = dlen + 1
[1843]1138       ENDDO rloop1
[1864]1139
[1839]1140       ALLOCATE( lrd(1:dlen), ard(1:dlen), crd(1:dlen) )
[1864]1141       DO jj = 1, dlen+1
1142          BACKSPACE ( 201, IOSTAT=ierrn )
1143       ENDDO             
1144       DO jj = 1, dlen
1145          READ ( 201, * ) lrd(jj), ard(jj), crd(jj) 
[1839]1146       ENDDO
[1819]1147
[1864]1148!
1149!--    Read tables (turb_cl(alpha),turb_cd(alpha) for the different profiles:
[1819]1150
[1864]1151       dlen = 0
1152
[1839]1153       READ ( 201, '(3/)' )
[1864]1154
[1843]1155       rloop2: DO
[1864]1156          READ ( 201, *, IOSTAT=ierrn ) chmess
1157          IF ( ierrn < 0  .OR.  chmess == '#'  .OR.  chmess == '')  EXIT rloop2
1158          dlen = dlen + 1
1159       ENDDO rloop2 
1160
[1912]1161       ALLOCATE( alpha_attack_tab(1:dlen), turb_cl_table(1:dlen,1:nairfoils),  &
[1864]1162                 turb_cd_table(1:dlen,1:nairfoils) )
1163
1164       DO jj = 1,dlen+1
1165          BACKSPACE ( 201, IOSTAT=ierrn )
1166       ENDDO 
1167
[1839]1168       DO jj = 1,dlen
[1864]1169          READ ( 201, * ) read_cl_cd
1170          alpha_attack_tab(jj) = read_cl_cd(1)
1171          DO ii= 1, nairfoils
1172             turb_cl_table(jj,ii) = read_cl_cd(ii*2)
1173             turb_cd_table(jj,ii) = read_cl_cd(ii*2+1)
[1819]1174          ENDDO
[1864]1175
[1839]1176       ENDDO
[1864]1177
1178       dlenbl = dlen 
1179
[1839]1180       CLOSE ( 201 )
[1819]1181
[1864]1182!
1183!--    For each possible radial position (resolution: 0.1 m --> 630 values) and
1184!--    each possible angle of attack (resolution: 0.01 degrees --> 36000 values!)
1185!--    determine the lift and drag coefficient by interpolating between the
1186!--    tabulated values of each table (interpolate to current angle of attack)
1187!--    and between the tables (interpolate to current radial position):
[1839]1188
[1864]1189       ALLOCATE( turb_cl_sel1(0:dlenbl) ) 
1190       ALLOCATE( turb_cl_sel2(0:dlenbl) ) 
1191       ALLOCATE( turb_cd_sel1(0:dlenbl) )
1192       ALLOCATE( turb_cd_sel2(0:dlenbl) )
[1819]1193
[1864]1194       radres     = INT( rr(1) * 10.0_wp ) + 1_iwp
1195       dlenbl_int = INT( 360.0_wp / accu_cl_cd_tab ) + 1_iwp
[1839]1196
[1864]1197
1198       ALLOCATE( turb_cl_tab(0:dlenbl_int,1:radres) )
1199       ALLOCATE( turb_cd_tab(0:dlenbl_int,1:radres) )
1200
1201
1202       DO iir = 1, radres ! loop over radius
1203
1204          DO iialpha = 1, dlenbl_int  ! loop over angles
1205
1206             cur_r = ( iir - 1_iwp ) * 0.1_wp             
1207             alpha_attack_i = -180.0_wp + REAL( iialpha-1 ) * accu_cl_cd_tab
1208             ialpha = 1
1209
[2152]1210             DO WHILE ( ( alpha_attack_i > alpha_attack_tab(ialpha) ) .AND. (ialpha <= dlen ) )
[1864]1211                ialpha = ialpha + 1
1212             ENDDO
[1819]1213!
[1864]1214!--          Find position in table
1215             lct = MINLOC( ABS( trad1 - cur_r ) )
1216!                lct(1) = lct(1)
[1819]1217
[1864]1218             IF ( ( trad1(lct(1)) - cur_r ) .GT. 0.0 )  THEN
1219                lct(1) = lct(1) - 1
1220             ENDIF
[1819]1221
[1864]1222             trow = lct(1)
1223!
1224!--          Calculate weights for interpolation
1225             weight_a = ( trad2(trow) - cur_r ) / ( trad2(trow) - trad1(trow) )
1226             weight_b = ( cur_r - trad1(trow) ) / ( trad2(trow) - trad1(trow) )
1227             t1 = ttoint1(trow)
1228             t2 = ttoint2(trow)
[1819]1229
[1912]1230             IF ( t1 .EQ. t2 ) THEN  ! if both are the same, the weights are NaN
1231                weight_a = 0.5_wp    ! then do interpolate in between same twice
1232                weight_b = 0.5_wp    ! using 0.5 as weight
[1864]1233             ENDIF
[1819]1234
[1864]1235             IF ( t1 == 0 .AND. t2 == 0 ) THEN
1236                turb_cd_sel1 = 0.0_wp
1237                turb_cd_sel2 = 0.0_wp
1238                turb_cl_sel1 = 0.0_wp
1239                turb_cl_sel2 = 0.0_wp
1240             ELSE
1241                turb_cd_sel1 = turb_cd_table(:,t1)
1242                turb_cd_sel2 = turb_cd_table(:,t2)
1243                turb_cl_sel1 = turb_cl_table(:,t1)
1244                turb_cl_sel2 = turb_cl_table(:,t2)
1245             ENDIF
1246
[1912]1247!
1248!--          Interpolation of lift and drag coefficiencts on fine grid of radius
1249!--          segments and angles of attack
1250
1251             turb_cl_tab(iialpha,iir) = ( alpha_attack_tab(ialpha) -           &
1252                                          alpha_attack_i ) /                   &
[1864]1253                                        ( alpha_attack_tab(ialpha) -           &
1254                                          alpha_attack_tab(ialpha-1) ) *       &
1255                                        ( weight_a * turb_cl_sel1(ialpha-1) +  &
1256                                          weight_b * turb_cl_sel2(ialpha-1) ) +&
1257                                        ( alpha_attack_i             -         &
1258                                          alpha_attack_tab(ialpha-1) ) /       &
1259                                        ( alpha_attack_tab(ialpha) -           &
1260                                          alpha_attack_tab(ialpha-1) ) *       &
1261                                        ( weight_a * turb_cl_sel1(ialpha) +    &
1262                                          weight_b * turb_cl_sel2(ialpha) )
[1912]1263             turb_cd_tab(iialpha,iir) = ( alpha_attack_tab(ialpha) -           &
1264                                          alpha_attack_i ) /                   &
[1864]1265                                        ( alpha_attack_tab(ialpha) -           &
1266                                          alpha_attack_tab(ialpha-1) ) *       &
1267                                        ( weight_a * turb_cd_sel1(ialpha-1) +  &
1268                                          weight_b * turb_cd_sel2(ialpha-1) ) +&
1269                                        ( alpha_attack_i             -         &
1270                                          alpha_attack_tab(ialpha-1) ) /       &
1271                                        ( alpha_attack_tab(ialpha) -           &
1272                                          alpha_attack_tab(ialpha-1) ) *       &
1273                                        ( weight_a * turb_cd_sel1(ialpha) +    &
1274                                          weight_b * turb_cd_sel2(ialpha) )
[1819]1275   
[1912]1276          ENDDO   ! end loop over angles of attack
[1819]1277
[1864]1278       ENDDO   ! end loop over radius
1279   
[1912]1280    END SUBROUTINE wtm_read_blade_tables
[1819]1281
1282
[1864]1283!------------------------------------------------------------------------------!
1284! Description:
1285! ------------
1286!> The projection matrix for the coordinate system of therotor disc in respect
1287!> to the yaw and tilt angle of the rotor is calculated
1288!------------------------------------------------------------------------------!
1289    SUBROUTINE wtm_rotate_rotor( inot )
[1819]1290
[1864]1291
1292       IMPLICIT NONE
1293
1294       INTEGER(iwp) :: inot
1295!
1296!--    Calculation of the rotation matrix for the application of the tilt to
1297!--    the rotors
1298       rot_eigen_rad(1) = SIN( phi_yaw(inot) )    ! x-component of the radial eigenvector
1299       rot_eigen_rad(2) = COS( phi_yaw(inot) )    ! y-component of the radial eigenvector
1300       rot_eigen_rad(3) = 0.0_wp                  ! z-component of the radial eigenvector
1301
1302       rot_eigen_azi(1) = 0.0_wp                  ! x-component of the azimuth eigenvector
1303       rot_eigen_azi(2) = 0.0_wp                  ! y-component of the azimuth eigenvector
1304       rot_eigen_azi(3) = 1.0_wp                  ! z-component of the azimuth eigenvector
1305
1306       rot_eigen_nor(1) =  COS( phi_yaw(inot) )   ! x-component of the normal eigenvector
1307       rot_eigen_nor(2) = -SIN( phi_yaw(inot) )   ! y-component of the normal eigenvector
1308       rot_eigen_nor(3) = 0.0_wp                  ! z-component of the normal eigenvector
[1839]1309   
[1864]1310!
1311!--    Calculation of the coordinate transformation matrix to apply a tilt to
1312!--    the rotor. If tilt = 0, rot_coord_trans is a unit matrix.
[1819]1313
[1912]1314       rot_coord_trans(inot,1,1) = rot_eigen_rad(1)**2                   *     &
[1864]1315                                   ( 1.0_wp - COS( tilt ) ) + COS( tilt ) 
[1912]1316       rot_coord_trans(inot,1,2) = rot_eigen_rad(1) * rot_eigen_rad(2)   *     &
1317                                   ( 1.0_wp - COS( tilt ) )              -     &
[1864]1318                                   rot_eigen_rad(3) * SIN( tilt )
[1912]1319       rot_coord_trans(inot,1,3) = rot_eigen_rad(1) * rot_eigen_rad(3)   *     &
1320                                   ( 1.0_wp - COS( tilt ) )              +     &
[1864]1321                                   rot_eigen_rad(2) * SIN( tilt )
[1912]1322       rot_coord_trans(inot,2,1) = rot_eigen_rad(2) * rot_eigen_rad(1)   *     &
1323                                   ( 1.0_wp - COS( tilt ) )              +     &
[1864]1324                                   rot_eigen_rad(3) * SIN( tilt )
[1912]1325       rot_coord_trans(inot,2,2) = rot_eigen_rad(2)**2                   *     &
[1864]1326                                   ( 1.0_wp - COS( tilt ) ) + COS( tilt ) 
[1912]1327       rot_coord_trans(inot,2,3) = rot_eigen_rad(2) * rot_eigen_rad(3)   *     &
1328                                   ( 1.0_wp - COS( tilt ) )              -     &
[1864]1329                                   rot_eigen_rad(1) * SIN( tilt )
[1912]1330       rot_coord_trans(inot,3,1) = rot_eigen_rad(3) * rot_eigen_rad(1)   *     &
1331                                   ( 1.0_wp - COS( tilt ) )              -     &
[1864]1332                                   rot_eigen_rad(2) * SIN( tilt )
[1912]1333       rot_coord_trans(inot,3,2) = rot_eigen_rad(3) * rot_eigen_rad(2)   *     &
1334                                   ( 1.0_wp - COS( tilt ) )              +     &
[1864]1335                                   rot_eigen_rad(1) * SIN( tilt )
[1912]1336       rot_coord_trans(inot,3,3) = rot_eigen_rad(3)**2                   *     &
[1864]1337                                   ( 1.0_wp - COS( tilt ) ) + COS( tilt )
[1839]1338
[1864]1339!
1340!--    Vectors for the Transformation of forces from the rotor's spheric
1341!--    coordinate system to the cartesian coordinate system
1342       rotx(inot,:) = MATMUL( rot_coord_trans(inot,:,:), rot_eigen_nor )
1343       roty(inot,:) = MATMUL( rot_coord_trans(inot,:,:), rot_eigen_rad )
1344       rotz(inot,:) = MATMUL( rot_coord_trans(inot,:,:), rot_eigen_azi )
1345   
1346    END SUBROUTINE wtm_rotate_rotor
[1839]1347
1348
[1819]1349!------------------------------------------------------------------------------!
1350! Description:
1351! ------------
[1839]1352!> Calculation of the forces generated by the wind turbine
[1819]1353!------------------------------------------------------------------------------!
1354    SUBROUTINE wtm_forces
1355
[1864]1356
[1819]1357       IMPLICIT NONE
1358
1359       CHARACTER (LEN=2) ::  turbine_id
1360
[1839]1361       INTEGER(iwp) ::  i, j, k          !< loop indices
1362       INTEGER(iwp) ::  inot             !< turbine loop index (turbine id)
1363       INTEGER(iwp) ::  iialpha, iir     !<
1364       INTEGER(iwp) ::  rseg, rseg_int   !<
1365       INTEGER(iwp) ::  ring, ring_int   !<
1366       INTEGER(iwp) ::  ii, jj, kk       !<
[1912]1367   
[2232]1368       REAL(wp)     ::  flag               !< flag to mask topography grid points
[1839]1369       REAL(wp)     ::  sin_rot, cos_rot   !<
1370       REAL(wp)     ::  sin_yaw, cos_yaw   !<
[1912]1371       
1372       REAL(wp) ::  aa, bb, cc, dd  !< interpolation distances
1373       REAL(wp) ::  gg              !< interpolation volume var 
1374       
1375       REAL(wp) ::  dist_u_3d, dist_v_3d, dist_w_3d  !< smearing distances
1376
1377       
[1839]1378!
[1819]1379!      Variables for pitch control
1380       REAL(wp)     ::  torque_max=0.0_wp
[1839]1381       LOGICAL      ::  pitch_sw=.FALSE.
1382
[1819]1383       INTEGER(iwp), DIMENSION(1) :: lct=0
[1839]1384       REAL(wp), DIMENSION(1)     :: rad_d=0.0_wp
[2152]1385       
1386       REAL(wp) :: tl_factor !< factor for tip loss correction
[1819]1387
1388
[1864]1389       CALL cpu_log( log_point_s(61), 'wtm_forces', 'start' )
[1819]1390
1391
1392       IF ( simulated_time >= time_turbine_on ) THEN
1393
[1864]1394!
[1819]1395!--       Set forces to zero for each new time step:
1396          thrust(:,:,:)         = 0.0_wp
1397          torque_y(:,:,:)       = 0.0_wp
1398          torque_z(:,:,:)       = 0.0_wp
1399          torque_total(:)       = 0.0_wp
[1912]1400          rot_tend_x(:,:,:)     = 0.0_wp
1401          rot_tend_y(:,:,:)     = 0.0_wp
1402          rot_tend_z(:,:,:)     = 0.0_wp
1403          thrust_rotor(:)       = 0.0_wp
[1819]1404!
1405!--       Loop over number of turbines:
1406          DO inot = 1, nturbines
1407
1408             cos_yaw = COS(phi_yaw(inot))
1409             sin_yaw = SIN(phi_yaw(inot))
1410!
[1839]1411!--          Loop over rings of each turbine:
[1819]1412             DO ring = 1, nrings(inot)
1413
1414                thrust_seg(:)   = 0.0_wp
1415                torque_seg_y(:) = 0.0_wp
1416                torque_seg_z(:) = 0.0_wp
1417!
1418!--             Determine distance between each ring (center) and the hub:
1419                cur_r = (ring - 0.5_wp) * delta_r(inot)
1420
1421!
[1839]1422!--             Loop over segments of each ring of each turbine:
[1819]1423                DO rseg = 1, nsegs(ring,inot)
1424!
[1864]1425!--                !-----------------------------------------------------------!
1426!--                !-- Determine coordinates of the ring segments            --!
1427!--                !-----------------------------------------------------------!
[1819]1428!
[1864]1429!--                Determine angle of ring segment towards zero degree angle of
1430!--                rotor system (at zero degree rotor direction vectors aligned
1431!--                with y-axis):
[1819]1432                   phi_rotor = rseg * 2.0_wp * pi / nsegs(ring,inot)
1433                   cos_rot   = COS( phi_rotor )
1434                   sin_rot   = SIN( phi_rotor )
1435
[1864]1436!--                Now the direction vectors can be determined with respect to
1437!--                the yaw and tilt angle:
[1819]1438                   re(1) = cos_rot * sin_yaw
[1839]1439                   re(2) = cos_rot * cos_yaw   
[1819]1440                   re(3) = sin_rot
1441
1442                   rote = MATMUL( rot_coord_trans(inot,:,:), re )
1443!
1444!--                Coordinates of the single segments (center points):
1445                   rbx(ring,rseg) = rcx(inot) + cur_r * rote(1)
1446                   rby(ring,rseg) = rcy(inot) + cur_r * rote(2)
1447                   rbz(ring,rseg) = rcz(inot) + cur_r * rote(3)
1448
[1864]1449!--                !-----------------------------------------------------------!
1450!--                !-- Interpolation of the velocity components from the     --!
1451!--                !-- cartesian grid point to the coordinates of each ring  --!
1452!--                !-- segment (follows a method used in the particle model) --!
1453!--                !-----------------------------------------------------------!
[1819]1454
1455                   u_int(inot,ring,rseg)     = 0.0_wp
1456                   u_int_1_l(inot,ring,rseg) = 0.0_wp
1457
1458                   v_int(inot,ring,rseg)     = 0.0_wp
1459                   v_int_1_l(inot,ring,rseg) = 0.0_wp
1460
1461                   w_int(inot,ring,rseg)     = 0.0_wp
1462                   w_int_1_l(inot,ring,rseg) = 0.0_wp
1463
1464!
1465!--                Interpolation of the u-component:
1466
1467                   ii =   rbx(ring,rseg) * ddx
1468                   jj = ( rby(ring,rseg) - 0.5_wp * dy ) * ddy
1469                   kk = ( rbz(ring,rseg) - 0.5_wp * dz ) / dz
1470!
[1864]1471!--                Interpolate only if all required information is available on
1472!--                the current PE:
[1839]1473                   IF ( ( ii >= nxl )  .AND.  ( ii <= nxr ) )  THEN
1474                      IF ( ( jj >= nys )  .AND.  ( jj <= nyn ) )  THEN
[1819]1475
[1839]1476                         aa = ( ( ii + 1          ) * dx - rbx(ring,rseg) ) *  &
1477                              ( ( jj + 1 + 0.5_wp ) * dy - rby(ring,rseg) )
1478                         bb = ( rbx(ring,rseg) - ii * dx )                  *  &
1479                              ( ( jj + 1 + 0.5_wp ) * dy - rby(ring,rseg) )
1480                         cc = ( ( ii+1            ) * dx - rbx(ring,rseg) ) *  &
1481                              ( rby(ring,rseg) - ( jj + 0.5_wp ) * dy )
1482                         dd = ( rbx(ring,rseg) -              ii * dx )     *  &
1483                              ( rby(ring,rseg) - ( jj + 0.5_wp ) * dy ) 
[1819]1484                         gg = dx * dy
1485
[1864]1486                         u_int_l = ( aa * u(kk,jj,ii)     +                    &
1487                                     bb * u(kk,jj,ii+1)   +                    &
1488                                     cc * u(kk,jj+1,ii)   +                    &
1489                                     dd * u(kk,jj+1,ii+1)                      &
[1819]1490                                   ) / gg
1491
[1864]1492                         u_int_u = ( aa * u(kk+1,jj,ii)     +                  &
1493                                     bb * u(kk+1,jj,ii+1)   +                  &
1494                                     cc * u(kk+1,jj+1,ii)   +                  &
1495                                     dd * u(kk+1,jj+1,ii+1)                    &
[1819]1496                                   ) / gg
1497
[1864]1498                         u_int_1_l(inot,ring,rseg) = u_int_l          +        &
1499                                     ( rbz(ring,rseg) - zu(kk) ) / dz *        &
[1819]1500                                     ( u_int_u - u_int_l )
1501
1502                      ELSE
1503                         u_int_1_l(inot,ring,rseg) = 0.0_wp
1504                      ENDIF
1505                   ELSE
1506                      u_int_1_l(inot,ring,rseg) = 0.0_wp
1507                   ENDIF
1508
1509
1510!
1511!--                Interpolation of the v-component:
1512                   ii = ( rbx(ring,rseg) - 0.5_wp * dx ) * ddx
[1839]1513                   jj =   rby(ring,rseg)                 * ddy
[1819]1514                   kk = ( rbz(ring,rseg) + 0.5_wp * dz ) / dz 
1515!
[1864]1516!--                Interpolate only if all required information is available on
1517!--                the current PE:
[1839]1518                   IF ( ( ii >= nxl )  .AND.  ( ii <= nxr ) )  THEN
1519                      IF ( ( jj >= nys )  .AND.  ( jj <= nyn ) )  THEN
[1819]1520
[1839]1521                         aa = ( ( ii + 1 + 0.5_wp ) * dx - rbx(ring,rseg) ) *  &
1522                              ( ( jj + 1 )          * dy - rby(ring,rseg) )
1523                         bb = ( rbx(ring,rseg)     - ( ii + 0.5_wp ) * dx ) *  &
1524                              ( ( jj + 1 ) * dy          - rby(ring,rseg) )
1525                         cc = ( ( ii + 1 + 0.5_wp ) * dx - rbx(ring,rseg) ) *  &
1526                              ( rby(ring,rseg)           -        jj * dy )
1527                         dd = ( rbx(ring,rseg)     - ( ii + 0.5_wp ) * dx ) *  &
1528                              ( rby(ring,rseg)           -        jj * dy )
[1819]1529                         gg = dx * dy
1530
[1864]1531                         v_int_l = ( aa * v(kk,jj,ii)     +                    &
1532                                     bb * v(kk,jj,ii+1)   +                    &
1533                                     cc * v(kk,jj+1,ii)   +                    &
1534                                     dd * v(kk,jj+1,ii+1)                      &
[1819]1535                                   ) / gg
1536
[1864]1537                         v_int_u = ( aa * v(kk+1,jj,ii)     +                  &
1538                                     bb * v(kk+1,jj,ii+1)   +                  &
1539                                     cc * v(kk+1,jj+1,ii)   +                  &
1540                                     dd * v(kk+1,jj+1,ii+1)                    &
[1819]1541                                  ) / gg
1542
[1864]1543                         v_int_1_l(inot,ring,rseg) = v_int_l +                 &
1544                                     ( rbz(ring,rseg) - zu(kk) ) / dz *        &
[1819]1545                                     ( v_int_u - v_int_l )
1546
1547                      ELSE
1548                         v_int_1_l(inot,ring,rseg) = 0.0_wp
1549                      ENDIF
1550                   ELSE
1551                      v_int_1_l(inot,ring,rseg) = 0.0_wp
1552                   ENDIF
1553
1554
1555!
1556!--                Interpolation of the w-component:
1557                   ii = ( rbx(ring,rseg) - 0.5_wp * dx ) * ddx
1558                   jj = ( rby(ring,rseg) - 0.5_wp * dy ) * ddy
[1839]1559                   kk =   rbz(ring,rseg)                 / dz
[1819]1560!
[1864]1561!--                Interpolate only if all required information is available on
1562!--                the current PE:
[1839]1563                   IF ( ( ii >= nxl )  .AND.  ( ii <= nxr ) )  THEN
1564                      IF ( ( jj >= nys )  .AND.  ( jj <= nyn ) )  THEN
[1819]1565
[1839]1566                         aa = ( ( ii + 1 + 0.5_wp ) * dx - rbx(ring,rseg) ) *  &
1567                              ( ( jj + 1 + 0.5_wp ) * dy - rby(ring,rseg) )
1568                         bb = ( rbx(ring,rseg)     - ( ii + 0.5_wp ) * dx ) *  &
1569                              ( ( jj + 1 + 0.5_wp ) * dy - rby(ring,rseg) )
1570                         cc = ( ( ii + 1 + 0.5_wp ) * dx - rbx(ring,rseg) ) *  &
1571                              ( rby(ring,rseg)     - ( jj + 0.5_wp ) * dy )
1572                         dd = ( rbx(ring,rseg)     - ( ii + 0.5_wp ) * dx ) *  &
1573                              ( rby(ring,rseg)     - ( jj + 0.5_wp ) * dy )
[1819]1574                         gg = dx * dy
1575
[1864]1576                         w_int_l = ( aa * w(kk,jj,ii)     +                    &
1577                                     bb * w(kk,jj,ii+1)   +                    &
1578                                     cc * w(kk,jj+1,ii)   +                    &
1579                                     dd * w(kk,jj+1,ii+1)                      &
[1819]1580                                   ) / gg
1581
[1864]1582                         w_int_u = ( aa * w(kk+1,jj,ii)     +                  &
1583                                     bb * w(kk+1,jj,ii+1)   +                  &
1584                                     cc * w(kk+1,jj+1,ii)   +                  &
1585                                     dd * w(kk+1,jj+1,ii+1)                    &
[1819]1586                                    ) / gg
1587
[1864]1588                         w_int_1_l(inot,ring,rseg) = w_int_l +                 &
1589                                     ( rbz(ring,rseg) - zw(kk) ) / dz *        &
[1819]1590                                     ( w_int_u - w_int_l )
1591                      ELSE
1592                         w_int_1_l(inot,ring,rseg) = 0.0_wp
1593                      ENDIF
1594                   ELSE
1595                      w_int_1_l(inot,ring,rseg) = 0.0_wp
1596                   ENDIF
1597
1598                ENDDO
1599             ENDDO
1600
1601          ENDDO
1602
1603!
1604!--       Exchange between PEs (information required on each PE):
1605#if defined( __parallel )
[1839]1606          CALL MPI_ALLREDUCE( u_int_1_l, u_int, nturbines * MAXVAL(nrings) *   &
1607                              MAXVAL(nsegs), MPI_REAL, MPI_SUM, comm2d, ierr )
1608          CALL MPI_ALLREDUCE( v_int_1_l, v_int, nturbines * MAXVAL(nrings) *   &
1609                              MAXVAL(nsegs), MPI_REAL, MPI_SUM, comm2d, ierr )
1610          CALL MPI_ALLREDUCE( w_int_1_l, w_int, nturbines * MAXVAL(nrings) *   &
1611                              MAXVAL(nsegs), MPI_REAL, MPI_SUM, comm2d, ierr )
[1819]1612#else
1613          u_int = u_int_1_l
1614          v_int = v_int_1_l
1615          w_int = w_int_1_l
1616#endif
1617
1618
1619!
1620!--       Loop over number of turbines:
[1912]1621
[1819]1622          DO inot = 1, nturbines
[1912]1623pit_loop: DO
[1819]1624
[1912]1625             IF ( pitch_sw )  THEN
[1839]1626                torque_total(inot) = 0.0_wp
[1912]1627                thrust_rotor(inot) = 0.0_wp
1628                pitch_add(inot)    = pitch_add(inot) + 0.25_wp
1629!                 IF ( myid == 0 ) PRINT*, 'Pitch', inot, pitch_add(inot)
1630             ELSE
1631                cos_yaw = COS(phi_yaw(inot))
1632                sin_yaw = SIN(phi_yaw(inot))
1633                IF ( pitch_control )  THEN
1634                   pitch_add(inot) = MAX(pitch_add_old(inot) - pitch_rate *    &
1635                                         dt_3d , 0.0_wp )
1636                ENDIF
[1819]1637             ENDIF
1638
[1839]1639!
1640!--          Loop over rings of each turbine:
[1819]1641             DO ring = 1, nrings(inot)
1642!
1643!--             Determine distance between each ring (center) and the hub:
1644                cur_r = (ring - 0.5_wp) * delta_r(inot)
[1839]1645!
1646!--             Loop over segments of each ring of each turbine:
[1819]1647                DO rseg = 1, nsegs(ring,inot)
1648!
[1839]1649!--                Determine angle of ring segment towards zero degree angle of
1650!--                rotor system (at zero degree rotor direction vectors aligned
1651!--                with y-axis):
[1819]1652                   phi_rotor = rseg * 2.0_wp * pi / nsegs(ring,inot)
1653                   cos_rot   = COS(phi_rotor)
1654                   sin_rot   = SIN(phi_rotor)
1655!
[1839]1656!--                Now the direction vectors can be determined with respect to
1657!--                the yaw and tilt angle:
[1819]1658                   re(1) = cos_rot * sin_yaw
1659                   re(2) = cos_rot * cos_yaw
1660                   re(3) = sin_rot
1661
[1864]1662!                  The current unit vector in azimuthal direction:                         
[1819]1663                   rea(1) = - sin_rot * sin_yaw
1664                   rea(2) = - sin_rot * cos_yaw
1665                   rea(3) =   cos_rot
1666
1667!
[1864]1668!--                To respect the yawing angle for the calculations of
1669!--                velocities and forces the unit vectors perpendicular to the
1670!--                rotor area in direction of the positive yaw angle are defined:
[1819]1671                   ren(1) =   cos_yaw
1672                   ren(2) = - sin_yaw
1673                   ren(3) = 0.0_wp
1674!
1675!--                Multiplication with the coordinate transformation matrix
[1864]1676!--                gives the final unit vector with consideration of the rotor
1677!--                tilt:
[1819]1678                   rote = MATMUL( rot_coord_trans(inot,:,:), re )
1679                   rota = MATMUL( rot_coord_trans(inot,:,:), rea )
1680                   rotn = MATMUL( rot_coord_trans(inot,:,:), ren )
1681!
1682!--                Coordinates of the single segments (center points):
1683                   rbx(ring,rseg) = rcx(inot) + cur_r * rote(1)
1684
1685                   rby(ring,rseg) = rcy(inot) + cur_r * rote(2)
1686
1687                   rbz(ring,rseg) = rcz(inot) + cur_r * rote(3)
1688
1689!
[1864]1690!--                !-----------------------------------------------------------!
1691!--                !-- Calculation of various angles and relative velocities --!
1692!--                !-----------------------------------------------------------!
[1819]1693!
[1864]1694!--                In the following the 3D-velocity field is projected its
[1819]1695!--                components perpedicular and parallel to the rotor area
1696!--                The calculation of forces will be done in the rotor-
1697!--                coordinates y' and z.
1698!--                The yaw angle will be reintroduced when the force is applied
1699!--                on the hydrodynamic equations
[1864]1700!
1701!--                Projection of the xy-velocities relative to the rotor area
1702!
[1819]1703!--                Velocity perpendicular to the rotor area:
[1864]1704                   u_rot = u_int(inot,ring,rseg)*rotn(1) +                     &
1705                   v_int(inot,ring,rseg)*rotn(2) +                             &
[1819]1706                   w_int(inot,ring,rseg)*rotn(3)
1707!
[1864]1708!--                Projection of the 3D-velocity vector in the azimuthal
1709!--                direction:
1710                   vtheta(rseg) = rota(1) * u_int(inot,ring,rseg) +            & 
1711                                  rota(2) * v_int(inot,ring,rseg) +            &
[1819]1712                                  rota(3) * w_int(inot,ring,rseg)
1713!
[1864]1714!--                Determination of the angle phi_rel between the rotor plane
1715!--                and the direction of the flow relative to the rotor:
[1819]1716
[1864]1717                   phi_rel(rseg) = ATAN( u_rot /                               &
1718                                         ( omega_rot(inot) * cur_r -           &
[1819]1719                                           vtheta(rseg) ) )
1720
1721!
[1864]1722!--                Interpolation of the local pitch angle from tabulated values
1723!--                to the current radial position:
[1819]1724
1725                   lct=minloc(ABS(cur_r-lrd))
1726                   rad_d=cur_r-lrd(lct)
1727                   
1728                   IF (cur_r == 0.0_wp) THEN
1729                      alpha_attack(rseg) = 0.0_wp
1730                   ELSE IF (cur_r >= lrd(size(ard))) THEN
[1864]1731                      alpha_attack(rseg) = ( ard(size(ard)) +                  &
1732                                             ard(size(ard)-1) ) / 2.0_wp
[1819]1733                   ELSE
[1864]1734                      alpha_attack(rseg) = ( ard(lct(1)) *  &
1735                                             ( ( lrd(lct(1)+1) - cur_r ) /     &
1736                                               ( lrd(lct(1)+1) - lrd(lct(1)) ) &
1737                                             ) ) + ( ard(lct(1)+1) *           &
1738                                             ( ( cur_r - lrd(lct(1)) ) /       &
1739                                               ( lrd(lct(1)+1) - lrd(lct(1)) ) ) )
[1819]1740                   ENDIF
1741
1742!
[1864]1743!--                In Fortran radian instead of degree is used as unit for all
1744!--                angles. Therefore, a transformation from angles given in
1745!--                degree to angles given in radian is necessary here:
1746                   alpha_attack(rseg) = alpha_attack(rseg) *                   &
[1819]1747                                        ( (2.0_wp*pi) / 360.0_wp )
1748!
[1864]1749!--                Substraction of the local pitch angle to obtain the local
1750!--                angle of attack:
[1819]1751                   alpha_attack(rseg) = phi_rel(rseg) - alpha_attack(rseg)
1752!
[1864]1753!--                Preliminary transformation back from angles given in radian
1754!--                to angles given in degree:
1755                   alpha_attack(rseg) = alpha_attack(rseg) *                   &
[1819]1756                                        ( 360.0_wp / (2.0_wp*pi) )
1757!
[1864]1758!--                Correct with collective pitch angle:
[2015]1759                   alpha_attack = alpha_attack - pitch_add(inot)
[1819]1760
1761!
[1864]1762!--                Determination of the magnitude of the flow velocity relative
1763!--                to the rotor:
[1912]1764                   vrel(rseg) = SQRT( u_rot**2 +                               &
1765                                      ( omega_rot(inot) * cur_r -              &
[1819]1766                                        vtheta(rseg) )**2 )
1767
1768!
[1864]1769!--                !-----------------------------------------------------------!
1770!--                !-- Interpolation of chord as well as lift and drag       --!
1771!--                !-- coefficients from tabulated values                    --!
1772!--                !-----------------------------------------------------------!
[1819]1773
1774!
[1864]1775!--                Interpolation of the chord_length from tabulated values to
1776!--                the current radial position:
[1819]1777
1778                   IF (cur_r == 0.0_wp) THEN
1779                      chord(rseg) = 0.0_wp
1780                   ELSE IF (cur_r >= lrd(size(crd))) THEN
[1864]1781                      chord(rseg) = (crd(size(crd)) + ard(size(crd)-1)) / 2.0_wp
[1819]1782                   ELSE
[1864]1783                      chord(rseg) = ( crd(lct(1)) *                            &
1784                            ( ( lrd(lct(1)+1) - cur_r ) /                      &
1785                              ( lrd(lct(1)+1) - lrd(lct(1)) ) ) ) +            &
1786                            ( crd(lct(1)+1) *                                  &
1787                            ( ( cur_r-lrd(lct(1)) ) /                          &
1788                              ( lrd(lct(1)+1) - lrd(lct(1)) ) ) )
[1819]1789                   ENDIF
1790
1791!
1792!--                Determine index of current angle of attack, needed for
[1864]1793!--                finding the appropriate interpolated values of the lift and
1794!--                drag coefficients (-180.0 degrees = 0, +180.0 degrees = 36000,
1795!--                so one index every 0.01 degrees):
1796                   iialpha = CEILING( ( alpha_attack(rseg) + 180.0_wp )        &
[1839]1797                                      * ( 1.0_wp / accu_cl_cd_tab ) )
[1819]1798!
1799!--                Determine index of current radial position, needed for
[1864]1800!--                finding the appropriate interpolated values of the lift and
1801!--                drag coefficients (one index every 0.1 m):
[1819]1802                   iir = CEILING( cur_r * 10.0_wp )
1803!
[1864]1804!--                Read in interpolated values of the lift and drag coefficients
1805!--                for the current radial position and angle of attack:
[1839]1806                   turb_cl(rseg) = turb_cl_tab(iialpha,iir)
1807                   turb_cd(rseg) = turb_cd_tab(iialpha,iir)
[1819]1808
1809!
[1864]1810!--                Final transformation back from angles given in degree to
1811!--                angles given in radian:
1812                   alpha_attack(rseg) = alpha_attack(rseg) *                   &
[1819]1813                                        ( (2.0_wp*pi) / 360.0_wp )
1814
[2152]1815                   IF ( tl_cor )  THEN
1816                   
1817!--                  Tip loss correction following Schito
1818!--                  Schito applies the tip loss correction only to the lift force
1819!--                  Therefore, the tip loss correction is only applied to the lift
1820!--                  coefficient and not to the drag coefficient in our case
1821!--                 
1822                     tl_factor = ( 2.0 / pi ) *                                &
1823                          ACOS( EXP( -1.0 * ( 3.0 * ( rr(inot) - cur_r ) /     &
1824                          ( 2.0 * cur_r * abs( sin( phi_rel(rseg) ) ) ) ) ) )
1825                         
1826                     turb_cl(rseg)  = tl_factor * turb_cl(rseg)                                 
1827                                 
1828                   END IF               
[1819]1829!
1830!--                !-----------------------------------------------------!
1831!--                !-- Calculation of the forces                       --!
1832!--                !-----------------------------------------------------!
1833
1834!
[1864]1835!--                Calculate the pre_factor for the thrust and torque forces:
[1819]1836
1837                   pre_factor = 0.5_wp * (vrel(rseg)**2) * 3.0_wp *  &
1838                                chord(rseg) * delta_r(inot) / nsegs(ring,inot)
1839
1840!
[1864]1841!--                Calculate the thrust force (x-component of the total force)
1842!--                for each ring segment:
1843                   thrust_seg(rseg) = pre_factor *                             &
1844                                      ( turb_cl(rseg) * COS(phi_rel(rseg)) +   &
1845                                        turb_cd(rseg) * SIN(phi_rel(rseg)) )
[1819]1846
1847!
[1864]1848!--                Determination of the second of the additional forces acting
1849!--                on the flow in the azimuthal direction: force vector as basis
1850!--                for torque (torque itself would be the vector product of the
1851!--                radius vector and the force vector):
1852                   torque_seg = pre_factor *                                   &
1853                                ( turb_cl(rseg) * SIN(phi_rel(rseg)) -         &
[1839]1854                                  turb_cd(rseg) * COS(phi_rel(rseg)) )
[1819]1855!
1856!--                Decomposition of the force vector into two parts:
[1864]1857!--                One acting along the y-direction and one acting along the
1858!--                z-direction of the rotor coordinate system:
[1819]1859
1860                   torque_seg_y(rseg) = -torque_seg * sin_rot
1861                   torque_seg_z(rseg) =  torque_seg * cos_rot
1862
[1912]1863!
1864!--                Add the segment thrust to the thrust of the whole rotor
1865                   thrust_rotor(inot) = thrust_rotor(inot) +                   &
1866                                        thrust_seg(rseg)                   
1867                   
1868
[1819]1869                   torque_total(inot) = torque_total(inot) + (torque_seg * cur_r)
1870
[1864]1871                ENDDO   !-- end of loop over ring segments
[1819]1872
1873!
[1864]1874!--             Restore the forces into arrays containing all the segments of
1875!--             each ring:
[1819]1876                thrust_ring(ring,:)   = thrust_seg(:)
1877                torque_ring_y(ring,:) = torque_seg_y(:)
1878                torque_ring_z(ring,:) = torque_seg_z(:)
1879
1880
[1864]1881             ENDDO   !-- end of loop over rings
[1819]1882
1883
[1864]1884             CALL cpu_log( log_point_s(62), 'wtm_controller', 'start' )
[1819]1885
[1912]1886             
1887             IF ( speed_control )  THEN
1888!
1889!--             Calculation of the current generator speed for rotor speed control
1890             
1891!                                     
1892!--             The acceleration of the rotor speed is calculated from
1893!--             the force balance of the accelerating torque
1894!--             and the torque of the rotating rotor and generator
1895                om_rate = ( torque_total(inot) * air_dens * gear_eff -         &
1896                            gear_ratio * torque_gen_old(inot) ) /              &
1897                          ( inertia_rot +                                      & 
1898                            gear_ratio * gear_ratio * inertia_gen ) * dt_3d
[1819]1899
[1912]1900!
1901!--             The generator speed is given by the product of gear gear_ratio
1902!--             and rotor speed
1903                omega_gen(inot) = gear_ratio * ( omega_rot(inot) + om_rate )     
1904             
1905             ENDIF
1906             
[1864]1907             IF ( pitch_control )  THEN
[1819]1908
[1912]1909!
1910!--             If the current generator speed is above rated, the pitch is not
1911!--             saturated and the change from the last time step is within the
1912!--             maximum pitch rate, then the pitch loop is repeated with a pitch
1913!--             gain
1914                IF ( (  omega_gen(inot)  > rated_genspeed   )  .AND.           &
1915                     ( pitch_add(inot) < 25.0_wp ) .AND.                       &
1916                     ( pitch_add(inot) < pitch_add_old(inot) +                 & 
1917                       pitch_rate * dt_3d  ) ) THEN
[1864]1918                   pitch_sw = .TRUE.
[1912]1919!
1920!--                Go back to beginning of pit_loop                   
1921                   CYCLE pit_loop
[1819]1922                ENDIF
[1912]1923               
1924!
1925!--             The current pitch is saved for the next time step
1926                pitch_add_old(inot) = pitch_add(inot)
[1864]1927                pitch_sw = .FALSE.
[1819]1928             ENDIF
[1912]1929             EXIT pit_loop             
1930          ENDDO pit_loop ! Recursive pitch control loop
[1819]1931
[1864]1932
[1819]1933!
[1864]1934!--          Call the rotor speed controller
1935             
[1819]1936             IF ( speed_control )  THEN
1937!
[1864]1938!--             Find processor at i_hub, j_hub             
[1912]1939                IF ( ( nxl <= i_hub(inot) )  .AND.  ( nxr >= i_hub(inot) ) )   &
1940                   THEN
1941                   IF ( ( nys <= j_hub(inot) )  .AND.  ( nyn >= j_hub(inot) ) )&
1942                      THEN
[1864]1943                      CALL wtm_speed_control( inot )
1944                   ENDIF
[1912]1945                ENDIF
[1864]1946                               
[1819]1947             ENDIF
1948
1949
[1864]1950             CALL cpu_log( log_point_s(62), 'wtm_controller', 'stop' )
[1819]1951
[1864]1952             CALL cpu_log( log_point_s(63), 'wtm_smearing', 'start' )
[1819]1953
1954
[1864]1955!--          !-----------------------------------------------------------------!
1956!--          !--                  Regularization kernel                      --!
1957!--          !-- Smearing of the forces and interpolation to cartesian grid  --!
1958!--          !-----------------------------------------------------------------!
[1819]1959!
[1864]1960!--          The aerodynamic blade forces need to be distributed smoothly on
1961!--          several mesh points in order to avoid singular behaviour
[1819]1962!
1963!--          Summation over sum of weighted forces. The weighting factor
[1864]1964!--          (calculated in user_init) includes information on the distance
1965!--          between the center of the grid cell and the rotor segment under
1966!--          consideration
[1819]1967!
[1864]1968!--          To save computing time, apply smearing only for the relevant part
1969!--          of the model domain:
[1819]1970!
1971!--
1972!--          Calculation of the boundaries:
[1864]1973             i_smear(inot) = CEILING( ( rr(inot) * ABS( roty(inot,1) ) +       &
1974                                        eps_min ) / dx )
1975             j_smear(inot) = CEILING( ( rr(inot) * ABS( roty(inot,2) ) +       &
1976                                        eps_min ) / dy )
[1819]1977
[1864]1978             DO i = MAX( nxl, i_hub(inot) - i_smear(inot) ),                   &
[1819]1979                    MIN( nxr, i_hub(inot) + i_smear(inot) )
[1864]1980                DO j = MAX( nys, j_hub(inot) - j_smear(inot) ),                &
[1819]1981                        MIN( nyn, j_hub(inot) + j_smear(inot) )
[2232]1982!                    DO k = MAX( nzb_u_inner(j,i)+1, k_hub(inot) - k_smear(inot) ), &
1983!                                 k_hub(inot) + k_smear(inot)
1984                   DO  k = nzb+1, k_hub(inot) + k_smear(inot)
[1819]1985                      DO ring = 1, nrings(inot)
1986                         DO rseg = 1, nsegs(ring,inot)
1987!
[2232]1988!--                         Predetermine flag to mask topography
1989                            flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 1 ) )
1990     
1991!
[1819]1992!--                         Determine the square of the distance between the
1993!--                         current grid point and each rotor area segment:
[1864]1994                            dist_u_3d = ( i * dx               - rbx(ring,rseg) )**2 + &
1995                                        ( j * dy + 0.5_wp * dy - rby(ring,rseg) )**2 + &
1996                                        ( k * dz - 0.5_wp * dz - rbz(ring,rseg) )**2
1997                            dist_v_3d = ( i * dx + 0.5_wp * dx - rbx(ring,rseg) )**2 + &
1998                                        ( j * dy               - rby(ring,rseg) )**2 + &
1999                                        ( k * dz - 0.5_wp * dz - rbz(ring,rseg) )**2
2000                            dist_w_3d = ( i * dx + 0.5_wp * dx - rbx(ring,rseg) )**2 + &
2001                                        ( j * dy + 0.5_wp * dy - rby(ring,rseg) )**2 + &
2002                                        ( k * dz               - rbz(ring,rseg) )**2
[1819]2003
2004!
2005!--                         3D-smearing of the forces with a polynomial function
2006!--                         (much faster than the old Gaussian function), using
2007!--                         some parameters that have been calculated in user_init.
2008!--                         The function is only similar to Gaussian function for
2009!--                         squared distances <= eps_min2:
2010                            IF ( dist_u_3d <= eps_min2 ) THEN
[2232]2011                               thrust(k,j,i) = thrust(k,j,i) +                     &
2012                                               thrust_ring(ring,rseg) *            &
2013                                               ( ( pol_a * dist_u_3d - pol_b ) *   & 
2014                                                dist_u_3d + 1.0_wp ) * eps_factor *&
2015                                                                       flag
[1819]2016                            ENDIF
2017                            IF ( dist_v_3d <= eps_min2 ) THEN
[2232]2018                               torque_y(k,j,i) = torque_y(k,j,i) +                    &
2019                                                 torque_ring_y(ring,rseg) *           &
2020                                                 ( ( pol_a * dist_v_3d - pol_b ) *    &
2021                                                  dist_v_3d + 1.0_wp ) * eps_factor * &
2022                                                                         flag
[1819]2023                            ENDIF
2024                            IF ( dist_w_3d <= eps_min2 ) THEN
[2232]2025                               torque_z(k,j,i) = torque_z(k,j,i) +                    &
2026                                                 torque_ring_z(ring,rseg) *           &
2027                                                 ( ( pol_a * dist_w_3d - pol_b ) *    &
2028                                                  dist_w_3d + 1.0_wp ) * eps_factor * &
2029                                                                         flag
[1819]2030                            ENDIF
2031
2032                         ENDDO  ! End of loop over rseg
2033                      ENDDO     ! End of loop over ring
2034             
2035!
2036!--                   Rotation of force components:
[2232]2037                      rot_tend_x(k,j,i) = rot_tend_x(k,j,i) + (                &
[1864]2038                                      thrust(k,j,i)*rotx(inot,1) +             &
2039                                      torque_y(k,j,i)*roty(inot,1) +           &
[2232]2040                                      torque_z(k,j,i)*rotz(inot,1)             &
2041                                                              ) * flag
[1819]2042                               
[2232]2043                      rot_tend_y(k,j,i) = rot_tend_y(k,j,i) + (                &
[1864]2044                                      thrust(k,j,i)*rotx(inot,2) +             &
2045                                      torque_y(k,j,i)*roty(inot,2) +           &
[2232]2046                                      torque_z(k,j,i)*rotz(inot,2)             &
2047                                                              ) * flag
[1819]2048                               
[2232]2049                      rot_tend_z(k,j,i) = rot_tend_z(k,j,i) + (                &
[1864]2050                                      thrust(k,j,i)*rotx(inot,3) +             &
2051                                      torque_y(k,j,i)*roty(inot,3) +           &
[2232]2052                                      torque_z(k,j,i)*rotz(inot,3)             &
2053                                                              ) * flag                   
[1819]2054
2055                   ENDDO        ! End of loop over k
2056                ENDDO           ! End of loop over j
2057             ENDDO              ! End of loop over i
2058
[1864]2059             CALL cpu_log( log_point_s(63), 'wtm_smearing', 'stop' )         
[1819]2060                   
[1912]2061          ENDDO                  !-- end of loop over turbines
[1819]2062
2063               
2064          IF ( yaw_control )  THEN
2065!
2066!--          Allocate arrays for yaw control at first call
2067!--          Can't be allocated before dt_3d is set
2068             IF ( start_up )  THEN
[1864]2069                WDLON = NINT( 30.0_wp / dt_3d )  ! 30s running mean array
[1819]2070                ALLOCATE( wd30(1:nturbines,1:WDLON) )
[1912]2071                wd30 = 999.0_wp                  ! Set to dummy value
[1819]2072                ALLOCATE( wd30_l(1:WDLON) )
2073               
[1864]2074                WDSHO = NINT( 2.0_wp / dt_3d )   ! 2s running mean array
[1819]2075                ALLOCATE( wd2(1:nturbines,1:WDSHO) )
[1912]2076                wd2 = 999.0_wp                   ! Set to dummy value
[1819]2077                ALLOCATE( wd2_l(1:WDSHO) )
2078                start_up = .FALSE.
2079             ENDIF         
2080
2081!
2082!--          Calculate the inflow wind speed
2083!--
2084!--          Loop over number of turbines:
2085             DO inot = 1, nturbines
2086!
2087!--             Find processor at i_hub, j_hub             
[1912]2088                IF ( ( nxl <= i_hub(inot) )  .AND.  ( nxr >= i_hub(inot) ) )   &
2089                   THEN
2090                   IF ( ( nys <= j_hub(inot) )  .AND.  ( nyn >= j_hub(inot) ) )&
2091                      THEN
[1864]2092
[1819]2093                      u_inflow_l(inot) = u(k_hub(inot),j_hub(inot),i_hub(inot))
[1864]2094
[1912]2095                      wdir_l(inot) = -1.0_wp * ATAN2(                          &
2096                         0.5_wp * ( v(k_hub(inot),j_hub(inot),i_hub(inot)+1) + &
2097                                    v(k_hub(inot),j_hub(inot),i_hub(inot)) ) , &
2098                         0.5_wp * ( u(k_hub(inot),j_hub(inot)+1,i_hub(inot)) + &
2099                                    u(k_hub(inot),j_hub(inot),i_hub(inot)) ) )
[1864]2100
2101                      CALL wtm_yawcontrol( inot )
2102
[1819]2103                      phi_yaw_l(inot) = phi_yaw(inot)
[2015]2104                     
[1819]2105                   ENDIF
2106                ENDIF
2107                   
[1864]2108             ENDDO                                 !-- end of loop over turbines
2109
2110!
[1929]2111!--          Transfer of information to the other cpus
2112#if defined( __parallel )         
[1864]2113             CALL MPI_ALLREDUCE( u_inflow_l, u_inflow, nturbines, MPI_REAL,    &
2114                                 MPI_SUM, comm2d, ierr )
2115             CALL MPI_ALLREDUCE( wdir_l, wdir, nturbines, MPI_REAL, MPI_SUM,   &
2116                                 comm2d, ierr )
2117             CALL MPI_ALLREDUCE( phi_yaw_l, phi_yaw, nturbines, MPI_REAL,      &
2118                                 MPI_SUM, comm2d, ierr )
[1929]2119#else
2120             u_inflow = u_inflow_l
2121             wdir     = wdir_l
2122             phi_yaw  = phi_yaw_l
[2015]2123             
2124             
[1929]2125#endif
[1819]2126             DO inot = 1, nturbines
[1864]2127!             
2128!--             Update rotor orientation               
2129                CALL wtm_rotate_rotor( inot )
[1819]2130
2131             ENDDO ! End of loop over turbines
2132                           
[2152]2133          END IF  ! end of yaw control
[1864]2134         
2135          IF ( speed_control )  THEN
2136!
2137!--          Transfer of information to the other cpus
[1912]2138!              CALL MPI_ALLREDUCE( omega_gen, omega_gen_old, nturbines,        &
2139!                                  MPI_REAL,MPI_SUM, comm2d, ierr )
[1929]2140#if defined( __parallel )   
[1864]2141             CALL MPI_ALLREDUCE( torque_gen, torque_gen_old, nturbines,        &
2142                                 MPI_REAL, MPI_SUM, comm2d, ierr )
2143             CALL MPI_ALLREDUCE( omega_rot_l, omega_rot, nturbines,            &
2144                                 MPI_REAL, MPI_SUM, comm2d, ierr )
[1912]2145             CALL MPI_ALLREDUCE( omega_gen_f, omega_gen_f_old, nturbines,      &
2146                                 MPI_REAL, MPI_SUM, comm2d, ierr )
[1929]2147#else
2148             torque_gen_old  = torque_gen
2149             omega_rot       = omega_rot_l
2150             omega_gen_f_old = omega_gen_f
2151#endif
[1864]2152           
2153          ENDIF
[1819]2154
2155          DO inot = 1, nturbines
2156
2157             IF ( myid == 0 ) THEN
2158                IF ( openfile_turb_mod(400+inot)%opened )  THEN
2159                   WRITE ( 400+inot, 106 ) simulated_time, omega_rot(inot),    &
[1912]2160                             omega_gen(inot), torque_gen_old(inot),            &
[1864]2161                             torque_total(inot), pitch_add(inot),              &
[1912]2162                             torque_gen_old(inot)*omega_gen(inot)*gen_eff,     &
[1864]2163                             torque_total(inot)*omega_rot(inot)*air_dens,      &
[1912]2164                             thrust_rotor(inot),                               & 
[1864]2165                             wdir(inot)*180.0_wp/pi,                           &
2166                             (phi_yaw(inot))*180.0_wp/pi                   
[1912]2167                             
[1819]2168                ELSE
2169
2170                   WRITE ( turbine_id,'(I2.2)')  inot
[1864]2171                   OPEN ( 400+inot, FILE=( 'TURBINE_PARAMETERS'//turbine_id ), &
2172                                            FORM='FORMATTED' )
[1819]2173                   WRITE ( 400+inot, 105 ) inot
2174                   WRITE ( 400+inot, 106 ) simulated_time, omega_rot(inot),    &
[1912]2175                             omega_gen(inot), torque_gen_old(inot),            &
[1864]2176                             torque_total(inot), pitch_add(inot),              &
[1912]2177                             torque_gen_old(inot)*omega_gen(inot)*gen_eff,     &
[1864]2178                             torque_total(inot)*omega_rot(inot)*air_dens,      &
[1912]2179                             thrust_rotor(inot),                               & 
[1864]2180                             wdir(inot)*180.0_wp/pi,                           &                   
2181                             (phi_yaw(inot))*180.0_wp/pi
[1819]2182                ENDIF
2183             ENDIF
2184
[1912]2185!--          Set open flag
[1819]2186             openfile_turb_mod(400+inot)%opened = .TRUE.
[1864]2187          ENDDO                                    !-- end of loop over turbines
[1819]2188
2189       ENDIF
2190
[1864]2191       CALL cpu_log( log_point_s(61), 'wtm_forces', 'stop' )
[1819]2192       
[1912]2193!
2194!--    Formats
[1819]2195       105 FORMAT ('Turbine control data for turbine ',I2,1X,':'/ &
2196              &'----------------------------------------'/ &
2197              &'   Time   RSpeed  GSpeed  ', &
[1912]2198               'GenTorque  AeroTorque  Pitch  Power(Gen)  Power(Rot)  ',       &
2199               'RotThrust  WDirection  YawOrient')
[1819]2200
[1864]2201       106 FORMAT (F9.3,2X,F7.3,2X,F7.2,2X,F9.1,3X,F9.1,1X,F6.2,2X,F10.1,2X,   &
2202                   F10.1,1X,F9.1,2X,F7.2,1X,F7.2)
[1819]2203
2204
2205    END SUBROUTINE wtm_forces
2206
2207   
[1839]2208!------------------------------------------------------------------------------!
2209! Description:
2210! ------------
2211!> Yaw controller for the wind turbine model
2212!------------------------------------------------------------------------------!
2213    SUBROUTINE wtm_yawcontrol( inot )
2214   
[1819]2215       USE constants
2216       USE kinds
[1839]2217               
2218       IMPLICIT NONE
[1819]2219     
2220       INTEGER(iwp)             :: inot
2221       INTEGER(iwp)             :: i_wd_30
2222       REAL(wp)                 :: missal
2223
2224       i_wd_30 = 0_iwp
2225
2226
[1864]2227!--    The yaw controller computes a 30s running mean of the wind direction.
2228!--    If the difference between turbine alignment and wind direction exceeds
2229!--    5°, the turbine is yawed. The mechanism stops as soon as the 2s-running
2230!--    mean of the missalignment is smaller than 0.5°.
2231!--    Attention: If the timestep during the simulation changes significantly
2232!--    the lengths of the running means change and it does not correspond to
2233!--    30s/2s anymore.
2234!--    ! Needs to be modified for these situations !
2235!--    For wind from the east, the averaging of the wind direction could cause
2236!--    problems and the yaw controller is probably flawed. -> Routine for
2237!--    averaging needs to be improved!
2238!
2239!--    Check if turbine is not yawing
[1819]2240       IF ( .NOT. doyaw(inot) )  THEN
[1843]2241!
[1864]2242!--       Write current wind direction into array
[1843]2243          wd30_l    = wd30(inot,:)
2244          wd30_l    = CSHIFT( wd30_l, SHIFT=-1 )
[1819]2245          wd30_l(1) = wdir(inot)
[1843]2246!
[1864]2247!--       Check if array is full ( no more dummies )
[1819]2248          IF ( .NOT. ANY( wd30_l == 999.) ) THEN
2249
2250             missal = SUM( wd30_l ) / SIZE( wd30_l ) - phi_yaw(inot)
2251!
[1864]2252!--          Check if turbine is missaligned by more than max_miss
[1843]2253             IF ( ABS( missal ) > max_miss )  THEN
2254!
[1864]2255!--             Check in which direction to yaw         
[1843]2256                yawdir(inot) = SIGN( 1.0_wp, missal )
[1819]2257!
[1864]2258!--             Start yawing of turbine
[1843]2259                phi_yaw(inot) = phi_yaw(inot) + yawdir(inot) * yaw_speed * dt_3d
[1819]2260                doyaw(inot) = .TRUE.
[1864]2261                wd30_l = 999.  ! fill with dummies again
[1819]2262             ENDIF
2263          ENDIF
2264         
2265          wd30(inot,:) = wd30_l
2266
2267!     
[1864]2268!--    If turbine is already yawing:
2269!--    Initialize 2 s running mean and yaw until the missalignment is smaller
2270!--    than min_miss
[1819]2271
2272       ELSE
2273!
2274!--       Initialize 2 s running mean
2275          wd2_l = wd2(inot,:)
[1864]2276          wd2_l = CSHIFT( wd2_l, SHIFT = -1 )
[1819]2277          wd2_l(1) = wdir(inot)
[1843]2278!     
[1864]2279!--       Check if array is full ( no more dummies )
2280          IF ( .NOT. ANY( wd2_l == 999.0_wp ) ) THEN
2281!
2282!--          Calculate missalignment of turbine       
[1819]2283             missal = SUM( wd2_l - phi_yaw(inot) ) / SIZE( wd2_l )
[1864]2284!
2285!--          Check if missalignment is still larger than 0.5 degree and if the
2286!--          yaw direction is still right
2287             IF ( ( ABS( missal ) > min_miss )  .AND.                          &
2288                  ( yawdir(inot) == SIGN( 1.0_wp, missal ) ) )  THEN
2289!
2290!--             Continue yawing       
2291                phi_yaw(inot) = phi_yaw(inot) + yawdir(inot) * yaw_speed * dt_3d
[1819]2292             ELSE
[1864]2293!
2294!--             Stop yawing       
[1819]2295                doyaw(inot) = .FALSE.
[1864]2296                wd2_l = 999.0_wp ! fill with dummies again
[1819]2297             ENDIF
2298          ELSE
[1864]2299!
2300!--          Continue yawing
[1843]2301             phi_yaw(inot) = phi_yaw(inot) + yawdir(inot) * yaw_speed * dt_3d
[1819]2302          ENDIF
2303     
2304          wd2(inot,:) = wd2_l
2305           
2306       ENDIF
2307     
[1839]2308    END SUBROUTINE wtm_yawcontrol 
[1819]2309
[1864]2310
[1819]2311!------------------------------------------------------------------------------!
2312! Description:
2313! ------------
[1864]2314!> Initialization of the speed control
2315!------------------------------------------------------------------------------!
2316    SUBROUTINE wtm_init_speed_control
2317
2318
2319       IMPLICIT NONE
2320
2321!
2322!--    If speed control is set, remaining variables and control_parameters for
2323!--    the control algorithm are calculated
2324!
2325!--    Calculate slope constant for region 15
2326       slope15   = ( slope2 * min_reg2 * min_reg2 ) / ( min_reg2 - min_reg15 )
2327!
2328!--    Calculate upper limit of slipage region
2329       vs_sysp   = rated_genspeed / 1.1_wp
2330!
2331!--    Calculate slope of slipage region
2332       slope25   = ( rated_power / rated_genspeed ) /                          &
2333                   ( rated_genspeed - vs_sysp )
2334!
2335!--    Calculate lower limit of slipage region
2336       min_reg25 = ( slope25 - SQRT( slope25 * ( slope25 - 4.0_wp *            &
2337                                                 slope2 * vs_sysp ) ) ) /      &
2338                   ( 2.0_wp * slope2 )
2339!
2340!--    Frequency for the simple low pass filter
2341       Fcorner   = 0.25_wp
2342!
2343!--    At the first timestep the torque is set to its maximum to prevent
2344!--    an overspeeding of the rotor
2345       torque_gen_old(:) = max_torque_gen 
2346     
2347    END SUBROUTINE wtm_init_speed_control
2348
2349
2350!------------------------------------------------------------------------------!
2351! Description:
2352! ------------
2353!> Simple controller for the regulation of the rotor speed
2354!------------------------------------------------------------------------------!
2355    SUBROUTINE wtm_speed_control( inot )
2356
2357
2358       IMPLICIT NONE
2359
[1912]2360       INTEGER(iwp)             :: inot
2361       
2362         
[1864]2363
2364!
2365!--    The controller is based on the fortran script from Jonkman
2366!--    et al. 2009 "Definition of a 5 MW Reference Wind Turbine for
2367!--    offshore system developement"
2368
2369!
2370!--    The generator speed is filtered by a low pass filter
2371!--    for the control of the generator torque       
2372       lp_coeff = EXP( -2.0_wp * 3.14_wp * dt_3d * Fcorner )
[1912]2373       omega_gen_f(inot) = ( 1.0_wp - lp_coeff ) * omega_gen(inot) + lp_coeff *&
[1864]2374                           omega_gen_f_old(inot)
2375
2376       IF ( omega_gen_f(inot) <= min_reg15 )  THEN
2377!                       
2378!--       Region 1: Generator torque is set to zero to accelerate the rotor:
2379          torque_gen(inot) = 0
2380       
2381       ELSEIF ( omega_gen_f(inot) <= min_reg2 )  THEN
2382!                       
2383!--       Region 1.5: Generator torque is increasing linearly with rotor speed:
2384          torque_gen(inot) = slope15 * ( omega_gen_f(inot) - min_reg15 )
2385                         
2386       ELSEIF ( omega_gen_f(inot) <= min_reg25 )  THEN
2387!
2388!--       Region 2: Generator torque is increased by the square of the generator
2389!--                 speed to keep the TSR optimal:
2390          torque_gen(inot) = slope2 * omega_gen_f(inot) * omega_gen_f(inot)
2391       
2392       ELSEIF ( omega_gen_f(inot) < rated_genspeed )  THEN
2393!                       
2394!--       Region 2.5: Slipage region between 2 and 3:
2395          torque_gen(inot) = slope25 * ( omega_gen_f(inot) - vs_sysp )
2396       
2397       ELSE
2398!                       
2399!--       Region 3: Generator torque is antiproportional to the rotor speed to
2400!--                 keep the power constant:
2401          torque_gen(inot) = rated_power / omega_gen_f(inot)
2402       
2403       ENDIF
2404!                       
2405!--    Calculate torque rate and confine with a max
2406       trq_rate = ( torque_gen(inot) - torque_gen_old(inot) ) / dt_3d
2407       trq_rate = MIN( MAX( trq_rate, -1.0_wp * max_trq_rate ), max_trq_rate )
2408!                       
2409!--    Calculate new gen torque and confine with max torque                         
2410       torque_gen(inot) = torque_gen_old(inot) + trq_rate * dt_3d
2411       torque_gen(inot) = MIN( torque_gen(inot), max_torque_gen )                                             
2412!
2413!--    Overwrite values for next timestep                       
[1912]2414       omega_rot_l(inot) = omega_gen(inot) / gear_ratio
[1864]2415
2416   
2417    END SUBROUTINE wtm_speed_control   
2418
2419
2420!------------------------------------------------------------------------------!
2421! Description:
2422! ------------
[1839]2423!> Application of the additional forces generated by the wind turbine on the
2424!> flow components (tendency terms)
2425!> Call for all grid points
[1819]2426!------------------------------------------------------------------------------!
2427    SUBROUTINE wtm_tendencies( component )
2428
[1839]2429   
[1819]2430       IMPLICIT NONE
2431
[1839]2432       INTEGER(iwp) ::  component   !< prognostic variable (u,v,w)
2433       INTEGER(iwp) ::  i           !< running index
2434       INTEGER(iwp) ::  j           !< running index
2435       INTEGER(iwp) ::  k           !< running index
[1819]2436
2437
2438       SELECT CASE ( component )
2439
2440       CASE ( 1 )
2441!
2442!--       Apply the x-component of the force to the u-component of the flow:
[1864]2443          IF ( simulated_time >= time_turbine_on )  THEN
[1819]2444             DO  i = nxlg, nxrg
2445                DO  j = nysg, nyng
[2232]2446                   DO  k = nzb+1, k_hub(1) + k_smear(1)
[1819]2447!
[1864]2448!--                   Calculate the thrust generated by the nacelle and the tower
[1912]2449                      tend_nac_x = 0.5_wp * nac_cd_surf(k,j,i) *               &
[1819]2450                                         SIGN( u(k,j,i)**2 , u(k,j,i) )     
[1912]2451                      tend_tow_x   = 0.5_wp * tow_cd_surf(k,j,i) *             &
[1819]2452                                         SIGN( u(k,j,i)**2 , u(k,j,i) ) 
2453                                                   
[2232]2454                      tend(k,j,i) = tend(k,j,i) + ( - rot_tend_x(k,j,i)        &
2455                                  - tend_nac_x - tend_tow_x )                  &
2456                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2457                                               BTEST( wall_flags_0(k,j,i), 1 ) )
[1819]2458                   ENDDO
2459                ENDDO
2460             ENDDO
2461          ENDIF
2462
2463       CASE ( 2 )
2464!
2465!--       Apply the y-component of the force to the v-component of the flow:
[1864]2466          IF ( simulated_time >= time_turbine_on )  THEN
[1819]2467             DO  i = nxlg, nxrg
2468                DO  j = nysg, nyng
[2232]2469                   DO  k = nzb+1, k_hub(1) + k_smear(1)
[1912]2470                      tend_nac_y = 0.5_wp * nac_cd_surf(k,j,i) *               &
[1819]2471                                         SIGN( v(k,j,i)**2 , v(k,j,i) )     
[1912]2472                      tend_tow_y   = 0.5_wp * tow_cd_surf(k,j,i) *             &
[1819]2473                                         SIGN( v(k,j,i)**2 , v(k,j,i) )                     
[2232]2474                      tend(k,j,i) = tend(k,j,i) + ( - rot_tend_y(k,j,i)        &
2475                                  - tend_nac_y - tend_tow_y )                  &
2476                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2477                                               BTEST( wall_flags_0(k,j,i), 2 ) )
[1819]2478                   ENDDO
2479                ENDDO
2480             ENDDO
2481          ENDIF
2482
2483       CASE ( 3 )
2484!
2485!--       Apply the z-component of the force to the w-component of the flow:
[1864]2486          IF ( simulated_time >= time_turbine_on )  THEN
[1819]2487             DO  i = nxlg, nxrg
2488                DO  j = nysg, nyng
[2232]2489                   DO  k = nzb+1,  k_hub(1) + k_smear(1)
2490                      tend(k,j,i) = tend(k,j,i) - rot_tend_z(k,j,i)            &
2491                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2492                                               BTEST( wall_flags_0(k,j,i), 3 ) )
[1819]2493                   ENDDO
2494                ENDDO
2495             ENDDO
2496          ENDIF
2497
2498
2499       CASE DEFAULT
2500
2501          WRITE( message_string, * ) 'unknown prognostic variable: ', component
2502          CALL message( 'wtm_tendencies', 'PA04??', 1, 2, 0, 6, 0 ) 
2503
2504       END SELECT
2505
2506
2507    END SUBROUTINE wtm_tendencies
2508
2509
2510!------------------------------------------------------------------------------!
2511! Description:
2512! ------------
[1839]2513!> Application of the additional forces generated by the wind turbine on the
2514!> flow components (tendency terms)
2515!> Call for grid point i,j
[1819]2516!------------------------------------------------------------------------------!
2517    SUBROUTINE wtm_tendencies_ij( i, j, component )
2518
2519
2520       IMPLICIT NONE
2521
[1839]2522       INTEGER(iwp) ::  component   !< prognostic variable (u,v,w)
2523       INTEGER(iwp) ::  i           !< running index
2524       INTEGER(iwp) ::  j           !< running index
2525       INTEGER(iwp) ::  k           !< running index
[1819]2526
2527       SELECT CASE ( component )
2528
2529       CASE ( 1 )
2530!
2531!--       Apply the x-component of the force to the u-component of the flow:
[1839]2532          IF ( simulated_time >= time_turbine_on )  THEN
[1819]2533
[2232]2534             DO  k = nzb+1,  k_hub(1) + k_smear(1)
[1819]2535!
[1839]2536!--             Calculate the thrust generated by the nacelle and the tower
[1912]2537                tend_nac_x = 0.5_wp * nac_cd_surf(k,j,i) *                     &
[1819]2538                                   SIGN( u(k,j,i)**2 , u(k,j,i) )     
[1912]2539                tend_tow_x   = 0.5_wp * tow_cd_surf(k,j,i) *                   &
[1819]2540                                   SIGN( u(k,j,i)**2 , u(k,j,i) ) 
[2232]2541                tend(k,j,i) = tend(k,j,i) + ( - rot_tend_x(k,j,i)              &
2542                            - tend_nac_x - tend_tow_x )                        &
2543                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2544                                               BTEST( wall_flags_0(k,j,i), 1 ) )
[1819]2545             ENDDO
2546          ENDIF
2547
2548       CASE ( 2 )
2549!
2550!--       Apply the y-component of the force to the v-component of the flow:
[1839]2551          IF ( simulated_time >= time_turbine_on )  THEN
[2232]2552             DO  k = nzb+1,  k_hub(1) + k_smear(1)
[1912]2553                tend_nac_y = 0.5_wp * nac_cd_surf(k,j,i) *                     &
[1819]2554                                   SIGN( v(k,j,i)**2 , v(k,j,i) )     
[1912]2555                tend_tow_y   = 0.5_wp * tow_cd_surf(k,j,i) *                   &
[1819]2556                                   SIGN( v(k,j,i)**2 , v(k,j,i) )                     
[2232]2557                tend(k,j,i) = tend(k,j,i) + ( - rot_tend_y(k,j,i)              &
2558                            - tend_nac_y - tend_tow_y )                        &
2559                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2560                                               BTEST( wall_flags_0(k,j,i), 2 ) )
[1819]2561             ENDDO
2562          ENDIF
2563
2564       CASE ( 3 )
2565!
2566!--       Apply the z-component of the force to the w-component of the flow:
[1839]2567          IF ( simulated_time >= time_turbine_on )  THEN
[2232]2568             DO  k = nzb+1,  k_hub(1) + k_smear(1)
2569                tend(k,j,i) = tend(k,j,i) - rot_tend_z(k,j,i)                  &
2570                                      * MERGE( 1.0_wp, 0.0_wp,                 &
2571                                               BTEST( wall_flags_0(k,j,i), 3 ) )
[1819]2572             ENDDO
2573          ENDIF
2574
2575
2576       CASE DEFAULT
2577
2578          WRITE( message_string, * ) 'unknown prognostic variable: ', component
2579          CALL message( 'wtm_tendencies', 'PA04??', 1, 2, 0, 6, 0 ) 
2580
2581       END SELECT
2582
2583
2584    END SUBROUTINE wtm_tendencies_ij
2585
2586 END MODULE wind_turbine_model_mod
Note: See TracBrowser for help on using the repository browser.