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

Last change on this file since 1929 was 1929, checked in by suehring, 8 years ago

several bugfixes in particle model and serial mode

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