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

Last change on this file since 2209 was 2153, checked in by lvollmer, 7 years ago

last commit documented

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