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

Last change on this file since 2101 was 2101, checked in by suehring, 7 years ago

last commit documented

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