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

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

Change unit number of file WTM_DATA

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