-
Notifications
You must be signed in to change notification settings - Fork 466
/
Copy pathMorison.f90
4318 lines (3635 loc) · 236 KB
/
Morison.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!**********************************************************************************************************************************
! The Morison and Morison_Types modules make up a template for creating user-defined calculations in the FAST Modularization
! Framework. Morison_Types will be auto-generated based on a description of the variables for the module.
!..................................................................................................................................
! LICENSING
! Copyright (C) 2012-2015 National Renewable Energy Laboratory
!
! This file is part of Morison.
!
! Licensed under the Apache License, Version 2.0 (the "License");
! you may not use this file except in compliance with the License.
! You may obtain a copy of the License at
!
! http://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
! See the License for the specific language governing permissions and
!
!**********************************************************************************************************************************
MODULE Morison
USE Waves
USE Morison_Types
USE Morison_Output
USE SeaSt_WaveField
! USE HydroDyn_Output_Types
USE NWTC_Library
USE YawOffset
IMPLICIT NONE
PRIVATE
TYPE(ProgDesc), PARAMETER :: Morison_ProgDesc = ProgDesc( 'Morison', '', '' )
! ..... Public Subroutines ...................................................................................................
PUBLIC:: Morison_GenerateSimulationNodes
PUBLIC :: Morison_Init ! Initialization routine
PUBLIC :: Morison_CalcOutput ! Routine for computing outputs
PUBLIC :: Morison_UpdateDiscState ! Routine for updating discrete states
CONTAINS
!----------------------------------------------------------------------------------------------------------------------------------
SUBROUTINE Morison_DirCosMtrx( pos0, pos1, DirCos )
! Compute the direction cosine matrix given two points along the axis of a cylinder
REAL(ReKi), INTENT( IN ) :: pos0(3), pos1(3)
Real(ReKi), INTENT( OUT ) :: DirCos(3,3)
Real(DbKi) :: xz, xyz
Real(DbKi) :: x0, y0, z0
Real(DbKi) :: x1, y1, z1
! Real(DbKi) :: temp
x0 = pos0(1)
y0 = pos0(2)
z0 = pos0(3)
x1 = pos1(1)
y1 = pos1(2)
z1 = pos1(3)
! Need to verify that z0 <= z1, but this was already handled in the element construction process!!! GJH 9/24/13
!IF ( z0 > z1 ) THEN
! temp = x0
! x0 = x1
! x1 = temp
! temp = y0
! y0 = y1
! y1 = temp
! temp = z0
! z0 = z1
! z1 = temp
!END IF
xz = sqrt((x0-x1)*(x0-x1)+(z0-z1)*(z0-z1))
xyz = sqrt((x0-x1)*(x0-x1)+(y0-y1)*(y0-y1)+(z0-z1)*(z0-z1))
IF ( xz==0 ) THEN
IF (y1<y0) THEN
DirCos = transpose(reshape((/ 1, 0, 0, 0, 0, -1, 0, 1, 0 /), shape(DirCos)))
ELSE
DirCos = transpose(reshape((/ 1, 0, 0, 0, 0, 1, 0, -1, 0 /), shape(DirCos)))
END IF
ELSE
DirCos(1, 1) = (z1-z0)/xz
DirCos(1, 2) = -(x1-x0)*(y1-y0)/(xz*xyz)
DirCos(1, 3) = (x1-x0)/xyz
DirCos(2, 1) = 0.0
DirCos(2, 2) = xz/xyz
DirCos(2, 3) = (y1-y0)/xyz
DirCos(3, 1) = -(x1-x0)/xz
DirCos(3, 2) = -(y1-y0)*(z1-z0)/(xz*xyz)
DirCos(3, 3) = (z1-z0)/xyz
END IF
END SUBROUTINE Morison_DirCosMtrx
!====================================================================================================
SUBROUTINE GetDistance ( a, b, l )
! This private subroutine computes the distance between points a and b.
!----------------------------------------------------------------------------------------------------
REAL(ReKi), INTENT ( IN ) :: a(3) ! the position of point a
REAL(ReKi), INTENT ( IN ) :: b(3) ! the position of point b
REAL(ReKi), INTENT ( OUT ) :: l ! the distance between point a and b
l = sqrt( ( a(1) - b(1) ) * ( a(1) - b(1) ) + ( a(2) - b(2) ) * ( a(2) - b(2) ) + ( a(3) - b(3) ) * ( a(3) - b(3) ) )
END SUBROUTINE GetDistance
!====================================================================================================
SUBROUTINE ElementCentroid ( Rs, Re, p1, h, DCM, centroid )
! This private subroutine computes the centroid of a tapered right cylinder element.
!----------------------------------------------------------------------------------------------------
REAL(ReKi), INTENT ( IN ) :: Rs ! starting radius
REAL(ReKi), INTENT ( IN ) :: Re ! ending radius
REAL(ReKi), INTENT ( IN ) :: p1(3) ! starting point of the element in global coordinates
REAL(ReKi), INTENT ( IN ) :: h ! height of the element
REAL(ReKi), INTENT ( IN ) :: DCM(3,3) ! direction cosine matrix to transform local element coordinates to global coordinates
REAL(ReKi), INTENT ( OUT ) :: centroid(3) ! centroid of the element in local coordinates
centroid(1) = 0.0
centroid(2) = 0.0
centroid(3) = h * (Rs*Rs + 2.0*Rs*Re + 3.0*Re*Re) / (4.0*( Rs*Rs + Rs*Re + Re*Re ) ) !( 2.0*Re + Rs ) / ( 3.0 * ( Rs + Re ) )
centroid = matmul( DCM, centroid ) + p1
END SUBROUTINE ElementCentroid
!====================================================================================================
REAL(ReKi) FUNCTION ElementVolume ( Rs, Re, h )
! This private function computes the volume of a tapered right cylinder element.
!----------------------------------------------------------------------------------------------------
REAL(ReKi), INTENT ( IN ) :: Rs ! starting radius
REAL(ReKi), INTENT ( IN ) :: Re ! ending radius
REAL(ReKi), INTENT ( IN ) :: h ! height of the element
ElementVolume = Pi*h*( Rs*Rs + Re*Re + Rs*Re ) / 3.0
END FUNCTION ElementVolume
!====================================================================================================
SUBROUTINE FindInterpFactor( p, p1, p2, s )
REAL(ReKi), INTENT ( IN ) :: p, p1, p2
REAL(ReKi), INTENT ( OUT ) :: s
REAL(ReKi) :: dp
! find normalized interpolation factor, s, such:
! p = p1*(1-s) + p2*s
! *--------------*--------------------------------*
! p1 p p2
!
! 0-----------------------------------------------1
! <------- s ---->
dp = p2 - p1
IF ( EqualRealNos(dp, 0.0_ReKi) ) THEN
s = 0
ELSE
s = ( p - p1 ) / dp
END IF
END SUBROUTINE FindInterpFactor
!=======================================================================
FUNCTION InterpWrappedStpInt( XValIn, XAry, YAry, Ind, AryLen )
! This funtion returns a y-value that corresponds to an input x-value which is wrapped back
! into the range [1-XAry(AryLen). It finds a x-value which corresponds to a value in the XAry where XAry(Ind-1) < MOD(XValIn, XAry(AryLen)) <= XAry(Ind)
! It is assumed that XAry is sorted in ascending order.
! It uses the passed index as the starting point and does a stepwise interpolation from there. This is
! especially useful when the calling routines save the value from the last time this routine was called
! for a given case where XVal does not change much from call to call. .
!
! This routine assumes YAry is INTEGER.
! Function declaration.
INTEGER :: InterpWrappedStpInt ! This function.
! Argument declarations.
INTEGER, INTENT(IN) :: AryLen ! Length of the arrays.
INTEGER, INTENT(INOUT) :: Ind ! Initial and final index into the arrays.
REAL(SiKi), INTENT(IN) :: XAry (AryLen) ! Array of X values to be interpolated.
REAL(SiKi), INTENT(IN) :: XValIn ! X value to be interpolated.
INTEGER, INTENT(IN) :: YAry (AryLen) ! Array of Y values to be interpolated.
REAL(SiKi) :: XVal ! X value to be interpolated.
! Wrap XValIn into the range XAry(1) to XAry(AryLen)
XVal = MOD(XValIn, XAry(AryLen))
! Set the Ind to the first index if we are at the beginning of XAry
IF ( XVal <= XAry(2) ) THEN
Ind = 1
END IF
! Let's check the limits first.
IF ( XVal <= XAry(1) ) THEN
InterpWrappedStpInt = YAry(1)
Ind = 1
RETURN
ELSE IF ( XVal >= XAry(AryLen) ) THEN
InterpWrappedStpInt = YAry(AryLen)
Ind = MAX(AryLen - 1, 1)
RETURN
END IF
! Let's interpolate!
Ind = MAX( MIN( Ind, AryLen-1 ), 1 )
DO
IF ( XVal < XAry(Ind) ) THEN
Ind = Ind - 1
ELSE IF ( XVal >= XAry(Ind+1) ) THEN
Ind = Ind + 1
ELSE
InterpWrappedStpInt = YAry(Ind)
RETURN
END IF
END DO
RETURN
END FUNCTION InterpWrappedStpInt ! ( XVal, XAry, YAry, Ind, AryLen )
!=======================================================================
FUNCTION InterpWrappedStpLogical( XValIn, XAry, YAry, Ind, AryLen )
! This funtion returns a y-value that corresponds to an input x-value which is wrapped back
! into the range [0-XAry(AryLen) by interpolating into the arrays.
! It is assumed that XAry is sorted in ascending order.
! It uses the passed index as the starting point and does a stepwise interpolation from there. This is
! especially useful when the calling routines save the value from the last time this routine was called
! for a given case where XVal does not change much from call to call. When there is no correlation
! from one interpolation to another, InterpBin() may be a better choice.
! It returns the first or last YAry() value if XVal is outside the limits of XAry().
! This routine assumes YAry is REAL.
! Function declaration.
LOGICAL :: InterpWrappedStpLogical ! This function.
! Argument declarations.
INTEGER, INTENT(IN) :: AryLen ! Length of the arrays.
INTEGER, INTENT(INOUT) :: Ind ! Initial and final index into the arrays.
REAL(SiKi), INTENT(IN) :: XAry (AryLen) ! Array of X values to be interpolated.
REAL(SiKi), INTENT(IN) :: XValIn ! X value to be interpolated.
LOGICAL, INTENT(IN) :: YAry (AryLen) ! Array of Y values to be interpolated.
REAL(SiKi) :: XVal ! X value to be interpolated.
! Wrap XValIn into the range XAry(1) to XAry(AryLen)
XVal = MOD(XValIn, XAry(AryLen))
! Set the Ind to the first index if we are at the beginning of XAry
IF ( XVal <= XAry(2) ) THEN
Ind = 1
END IF
! Let's check the limits first.
IF ( XVal <= XAry(1) ) THEN
InterpWrappedStpLogical = YAry(1)
Ind = 1
RETURN
ELSE IF ( XVal >= XAry(AryLen) ) THEN
InterpWrappedStpLogical = YAry(AryLen)
Ind = MAX(AryLen - 1, 1)
RETURN
END IF
! Let's interpolate!
Ind = MAX( MIN( Ind, AryLen-1 ), 1 )
DO
IF ( XVal < XAry(Ind) ) THEN
Ind = Ind - 1
ELSE IF ( XVal >= XAry(Ind+1) ) THEN
Ind = Ind + 1
ELSE
InterpWrappedStpLogical = YAry(Ind)
RETURN
END IF
END DO
RETURN
END FUNCTION InterpWrappedStpLogical ! ( XVal, XAry, YAry, Ind, AryLen )
!----------------------------------------------------------------------------------------------------------------------------------
subroutine GetOrientationAngles(p1, p2, phi, sinPhi, cosPhi, tanPhi, sinBeta, cosBeta, k_hat, errStat, errMsg)
real(ReKi), intent(in ) :: p1(3),p2(3)
real(ReKi), intent( out) :: phi, sinPhi, cosPhi, tanPhi, sinBeta, cosBeta, k_hat(3)
integer, intent( out) :: errStat ! returns a non-zero value when an error occurs
character(*), intent( out) :: errMsg ! Error message if errStat /= ErrID_None
character(*), parameter :: RoutineName = 'GetOrientationAngles'
real(ReKi) :: vec(3), vecLen, vecLen2D, beta
! Initialize errStat
errStat = ErrID_None
errMsg = ""
! calculate isntantaneous incline angle and heading, and related trig values
! the first and last NodeIndx values point to the corresponding Joint nodes idices which are at the start of the Mesh
vec = p2 - p1
vecLen = SQRT(Dot_Product(vec,vec))
vecLen2D = SQRT(vec(1)**2+vec(2)**2)
if ( vecLen < 0.000001 ) then
call SeterrStat(ErrID_Fatal, 'An element of the Morison structure has co-located endpoints! This should never occur. Please review your model.', errStat, errMsg, RoutineName )
return
else
k_hat = vec / vecLen
phi = atan2(vecLen2D, vec(3)) ! incline angle
end if
if ( EqualRealNos(phi, 0.0_ReKi) ) then
beta = 0.0_ReKi
else
beta = atan2(vec(2), vec(1)) ! heading of incline
endif
sinPhi = sin(phi)
cosPhi = cos(phi)
tanPhi = tan(phi)
sinBeta = sin(beta)
cosBeta = cos(beta)
end subroutine GetOrientationAngles
!----------------------------------------------------------------------------------------------------------------------------------
!function to return conical taper geometry calculations (volume and center of volume)
SUBROUTINE TaperCalc(R1, R2, H, taperV, h_c)
REAL(ReKi), INTENT ( IN ) :: R1
REAL(ReKi), INTENT ( IN ) :: R2
REAL(ReKi), INTENT ( IN ) :: H
REAL(ReKi), INTENT ( OUT ) :: taperV ! volume of tapered section
REAL(ReKi), INTENT ( OUT ) :: h_c ! center of mass offset from first node
real(ReKi) :: m
m = (R2-R1)/H
if ( EqualRealNos(R1, R2) ) then ! if just a cylinder
taperV = abs(pi*R1*R1*H)
h_c = H/2.0
elseif ( EqualRealNos(R1,0.0_ReKi) ) then ! seperate this case out because it gives a divide by zero in general formula
taperV = abs(1.0/3.0*pi*R2*R2*H) ! cone volume
h_c = 3.0/4.0*H ! from base
else
taperV = abs(pi/3.0/m*(R2**3 - R1**3))
h_c = H*(R1**2 + 2*R1*R2 + 3*R2**2)/4.0/(R1**2 + R1*R2 + R2**2) !( coneV*1./4.*coneH - coneVtip*(1./4.*(coneH-H) + H) )/ taperV ! from base
end if
END SUBROUTINE TaperCalc
!----------------------------------------------------------------------------------------------------------------------------------
SUBROUTINE CylInertia(R1, R2, H, rho, Il, Ir)
REAL(ReKi), INTENT ( IN ) :: R1
REAL(ReKi), INTENT ( IN ) :: R2
REAL(ReKi), INTENT ( IN ) :: H
REAL(ReKi), INTENT ( IN ) :: rho ! density of material
REAL(ReKi), INTENT ( OUT ) :: Il
REAL(ReKi), INTENT ( OUT ) :: Ir
real(ReKi) :: m, Ir_tip, h_c
m = (R2-R1)/H
if ( EqualRealNos(R1, R2) ) then ! if just a cylinder
Ir = abs(1.0/12.0* rho*pi*R1*R1*H *(3.0*R1*R1 + 4.0*H*H)) ! radial inertia about node 1
Il = abs(0.5* rho*pi*R1*R1*H *R1*R1)
ELSEIF ( EqualRealNos(R1,0.0_ReKi) ) then ! seperate this case out because it gives a divide by zero in general formula
Ir = abs(rho*pi*(1.0/20.0/m + 1.0/5.0/m**3) * R2**5)
Il = abs(1.0/10.0*rho*pi/m*R2**5)
ELSE
h_c = H*(R1**2 + 2*R1*R2 + 3*R2**2)/4.0/(R1**2 + R1*R2 + R2**2)
!l_c = R1/M + (R2-R1)/m *(R1**2 + 2*R1*R2 + 3*R2**2)/4/(R1**2 + R1*R2 + R2**2)
Ir_tip = abs(pi/20.0 *rho/m*(1.0 + 4.0/m**2) * (R2**5 - R1**5)) ! radial moment of inertia about tip of cone
Ir = abs(Ir_tip - rho/3.0/m*pi*(R2**3-R1**3) * (R1/m + 2.0*h_c)*R1/m ) ! radial moment of inertia about node 1
Il = abs(1.0/10.0/m*rho*pi*(R2**5 - R1**5))
END IF
END SUBROUTINE CylInertia
!----------------------------------------------------------------------------------------------------------------------------------
SUBROUTINE MarineGrowthPartSegment(R1, R2, Rmg1, Rmg2, L, rho, Vinner, Vouter, m_mg, h_c, Ilmg, Irmg)
REAL(ReKi), INTENT ( IN ) :: R1
REAL(ReKi), INTENT ( IN ) :: R2
REAL(ReKi), INTENT ( IN ) :: Rmg1
REAL(ReKi), INTENT ( IN ) :: Rmg2
REAL(ReKi), INTENT ( IN ) :: L
REAL(ReKi), INTENT ( IN ) :: rho ! density of material
REAL(ReKi), INTENT ( OUT ) :: Vinner ! volume from inner radius
REAL(ReKi), INTENT ( OUT ) :: Vouter ! volume from outer radius
REAL(ReKi), INTENT ( OUT ) :: m_mg ! mass of marine growth
REAL(ReKi), INTENT ( OUT ) :: h_c ! center of mass offset from first node
REAL(ReKi), INTENT ( OUT ) :: Ilmg ! moment of inertia about axis
REAL(ReKi), INTENT ( OUT ) :: Irmg ! moment of inertia about radial axis from first node
! Local variables
REAL(ReKi) :: cVinner ! center of volume from inner radius
REAL(ReKi) :: cVouter ! center of volume from outer radius
REAL(ReKi) :: Ilinner
REAL(ReKi) :: Irinner
REAL(ReKi) :: Ilouter
REAL(ReKi) :: Irouter
! get V and CV for element
call TaperCalc(R1, R2, L, Vinner, cVinner)
! get V and CV for marine growth displacement
call TaperCalc(Rmg1, Rmg2, L, Vouter, cVouter)
! get mass and CV specific to marine growth thickness
m_mg = (Vouter - Vinner)*rho
if ( EqualRealNos(m_mg, 0.0_ReKi) ) then
h_c = 0.0
else
h_c = (cVouter*Vouter - Vinner*cVinner)/(Vouter - Vinner)
end if
! get two moments of inertia for marine growth as if solid...
call CylInertia(Rmg1, Rmg2, L, rho, Ilouter, Irouter) ! inertias for marine growth if solid
call CylInertia(R1 , R2 , L, rho, Ilinner, Irinner) ! inertias for element if filled with marine growth
! subtract to get moments of inertia of marine growth shell
Ilmg = Ilouter - Ilinner
Irmg = Irouter - Irinner
END SUBROUTINE MarineGrowthPartSegment
!----------------------------------------------------------------------------------------------------------------------------------
SUBROUTINE FloodedBallastPartSegment(R1, R2, L, rho, V, m, h_c, Il, Ir)
REAL(ReKi), INTENT ( IN ) :: R1 ! interior radius of element at node point
REAL(ReKi), INTENT ( IN ) :: R2 ! interior radius of other end of part-element
REAL(ReKi), INTENT ( IN ) :: L ! distance positive along axis to end of part-element
REAL(ReKi), INTENT ( OUT ) :: V ! volume from inner radius
REAL(ReKi), INTENT ( IN ) :: rho ! density of ballast
REAL(ReKi), INTENT ( OUT ) :: m ! mass of material
REAL(ReKi), INTENT ( OUT ) :: h_c ! center of mass offset from first node
REAL(ReKi), INTENT ( OUT ) :: Il ! moment of inertia about axis
REAL(ReKi), INTENT ( OUT ) :: Ir ! moment of inertia about radial axis from first node
! get V and CV for flooded part of part-element
call TaperCalc(R1, R2, L, V, h_c)
m = rho*V
call CylInertia(R1, R2, L, rho, Il, Ir) ! inertias for filled section
END SUBROUTINE FloodedBallastPartSegment
!----------------------------------------------------------------------------------------------------------------------------------
SUBROUTINE WriteSummaryFile( UnSum, numJoints, numNodes, nodes, numMembers, members, &
NOutputs, OutParam, MOutLst, JOutLst, uMesh, yMesh, p, m, errStat, errMsg )
INTEGER, INTENT ( IN ) :: UnSum
INTEGER, INTENT ( IN ) :: numJoints
INTEGER, INTENT ( IN ) :: numNodes
TYPE(Morison_NodeType), ALLOCATABLE, INTENT ( IN ) :: nodes(:)
INTEGER, INTENT ( IN ) :: numMembers
TYPE(Morison_MemberType), ALLOCATABLE, INTENT ( IN ) :: members(:)
INTEGER, INTENT ( IN ) :: NOutputs
TYPE(OutParmType), ALLOCATABLE, INTENT ( IN ) :: OutParam(:)
TYPE(Morison_MOutput), ALLOCATABLE, INTENT ( IN ) :: MOutLst(:)
TYPE(Morison_JOutput), ALLOCATABLE, INTENT ( IN ) :: JOutLst(:)
TYPE(MeshType), INTENT ( INOUT ) :: uMesh
TYPE(MeshType), INTENT ( INOUT ) :: yMesh
TYPE(Morison_ParameterType), INTENT ( IN ) :: p
TYPE(Morison_MiscVarType), INTENT ( IN ) :: m
INTEGER, INTENT ( OUT ) :: errStat ! returns a non-zero value when an error occurs
CHARACTER(*), INTENT ( OUT ) :: errMsg ! Error message if errStat /= ErrID_None
INTEGER :: I, J, II
LOGICAL :: filledFlag ! flag indicating if element is filled/flooded
REAL(ReKi) :: ident(3,3) ! identity matrix
REAL(ReKi) :: ExtBuoyancy(6) ! sum of all external buoyancy forces lumped at (0,0,0)
REAL(ReKi) :: IntBuoyancy(6) ! sum of all internal buoyancy forces lumped at (0,0,0)
REAL(ReKi) :: MG_Wt(6) ! weight of the marine growth as applied to (0,0,0)
TYPE(MeshType) :: WRP_Mesh ! mesh representing the WAMIT reference point (0,0,0)
TYPE(MeshType) :: WRP_Mesh_position ! mesh representing the WAMIT reference point (0,0,0) (with no displaced position)
TYPE(MeshMapType) :: M_P_2_P ! Map Morison Point to WRP_Mesh point
REAL(ReKi) :: totalDisplVol ! total displaced volume of the structure
REAL(ReKi) :: totalVol ! total volume of structure
REAL(ReKi) :: MGvolume ! volume of the marine growth material
REAL(ReKi) :: totalMGVol !
REAL(ReKi) :: totalFillVol !
REAL(ReKi) :: COB(3) ! center of buoyancy location in global coordinates
INTEGER :: m1 ! Indices of the markers which surround the requested output location
REAL(ReKi) :: s ! The linear interpolation factor for the requested location
REAL(ReKi) :: outloc(3) ! Position of the requested member output
real(ReKi) :: pos(3), pos2(3) ! Position of a node or joint in the MSL inertial system
INTEGER :: mbrIndx, nodeIndx, c, N
CHARACTER(ChanLen) :: tmpName
REAL(ReKi) :: totalFillMass, mass_fill, memberVol
REAL(ReKi) :: totalMGMass
TYPE(Morison_NodeType) :: node1, node2
real(ReKi) :: ptLoad(6)
logical :: fillFlag
type(Morison_MemberType) :: mem
REAL(ReKi) :: Cd1, Cd2, Ca1, Ca2, Cp1, Cp2, AxCd1, AxCd2, AxCa1, AxCa2, AxCp1, AxCp2, Cb1, Cb2, JAxCd1, JAxCd2, JAxCa1, JAxCa2, JAxCp1, JAxCp2 ! tmp coefs
real(ReKi) :: F_B(6, numNodes), F_BF(6, numNodes), F_WMG(6, numNodes)
INTEGER :: ErrStat2
CHARACTER(ErrMsgLen) :: ErrMsg2
CHARACTER(*), PARAMETER :: RoutineName = 'WriteSummaryFile'
! Initialize data
errStat = ErrID_None
errMsg = ""
IF ( UnSum <= 0 ) RETURN ! can't write to the file (no summary file requested)
ExtBuoyancy = 0.0
totalFillMass = 0.0
totalDisplVol = 0.0
totalVol = 0.0
totalMGVol = 0.0
totalFillVol = 0.0
totalMGMass = 0.0
COB = 0.0
F_B = 0.0
F_BF = 0.0
F_WMG = 0.0
! Create identity matrix
CALL EYE(ident,errStat,errMsg)
do j = 1, numMembers
mem = members(j)
totalVol = totalVol + mem%Vouter
totalMGVol = totalMGVol + mem%Vouter - mem%Vinner
totalDisplVol = totalDisplVol + mem%Vsubmerged
totalFillVol = totalFillVol + mem%Vballast
do i = 1, mem%NElements
totalMGMass = totalMGMass + mem%m_mg_l(i)
totalMGMass = totalMGMass + mem%m_mg_u(i)
end do
do i = 1, mem%NElements+1
F_B (:,mem%NodeIndx(i)) = F_B (:,mem%NodeIndx(i)) + m%memberLoads(j)%F_B (:,i)
F_BF (:,mem%NodeIndx(i)) = F_BF (:,mem%NodeIndx(i)) + m%memberLoads(j)%F_BF (:,i)
F_WMG(:,mem%NodeIndx(i)) = F_WMG(:,mem%NodeIndx(i)) + m%memberLoads(j)%F_WMG(:,i)
end do
end do
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A37)' ) 'Strip-Theory Volume Calculations(m^3)'
WRITE( UnSum, '(A37)' ) '-------------------------------------'
WRITE( UnSum, '(A27,ES12.5)' ) ' Structure Volume : ', totalVol
WRITE( UnSum, '(A27,ES12.5)' ) ' Submerged Volume : ', totalDisplVol
WRITE( UnSum, '(A27,ES12.5)' ) ' Marine Growth Volume : ', totalMGVol
WRITE( UnSum, '(A27,ES12.5)' ) ' Ballasted Volume : ', totalFillVol
WRITE( UnSum, '(A111)') ' NOTE: Structure, Submerged and Marine Growth volumes are based on members not modelled with WAMIT'
WRITE( UnSum, '(A149)') ' Ballasted volume is computed from all members which are marked as filled in the HydroDyn input file, regardless of PropPot flag'
! Create a point mesh at (0,0,0) so that we can integrate the Morison load contributions to a single point for reporting purposes
CALL MeshCreate( BlankMesh = WRP_Mesh &
,IOS = COMPONENT_INPUT &
,Nnodes = 1 &
,errStat = errStat2 &
,ErrMess = errMsg2 &
,Force = .TRUE. &
,Moment = .TRUE. &
)
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
! Create the node on the mesh
CALL MeshPositionNode (WRP_Mesh &
, 1 &
, (/0.0_ReKi, 0.0_ReKi, 0.0_ReKi/) &
, errStat2 &
, errMsg2 &
)
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
IF ( errStat >= AbortErrLev ) then
call cleanup()
RETURN
end if
! Create the mesh element
CALL MeshConstructElement ( WRP_Mesh &
, ELEMENT_POINT &
, errStat2 &
, errMsg2 &
, 1 &
)
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
CALL MeshCommit ( WRP_Mesh &
, errStat2 &
, errMsg2 )
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
! we need the translation displacement mesh for loads transfer:
CALL MeshCopy ( SrcMesh = WRP_Mesh &
, DestMesh = WRP_Mesh_position &
, CtrlCode = MESH_SIBLING &
, IOS = COMPONENT_INPUT &
, TranslationDisp = .TRUE. &
, errStat = errStat2 &
, ErrMess = errMsg2 ) ! automatically sets DestMesh%RemapFlag = .TRUE.
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
IF ( errStat >= AbortErrLev ) then
call cleanup()
RETURN
end if
WRP_Mesh_position%TranslationDisp = 0.0 ! bjj: this is actually initialized in the ModMesh module, but I'll do it here anyway.
WRP_Mesh%RemapFlag = .TRUE.
! Attach the external distributed buoyancy loads to the distributed mesh so they can be transferred to the WRP
! Because of wave stretching and user-supplied waves, we may have loads above the still water line (SWL) which will be used
! in the hydrodynamics for conditions where the wave height is > SWL. So we now need to check that the vertical position
! is <= SWL for this summary file calculation.
DO J = 1, yMesh%Nnodes
if ( yMesh%Position(3,J) <= p%WaveField%MSL2SWL ) then ! need to check relative to MSL2SWL offset because the Mesh Positons are relative to MSL
if (J <= numJoints) then
ptLoad = F_B(:,J) + m%F_B_end(:,J)
else
ptLoad = F_B(:,J)
end if
yMesh%Force(:,J) = ptLoad(1:3)
yMesh%Moment(:,J) = ptLoad(4:6)
else
yMesh%Force(:,J) = 0.0
yMesh%Moment(:,J) = 0.0
end if ! <= still water line check
END DO ! DO J
! Transfer the loads from the distributed mesh to the (0,0,0) point mesh
CALL MeshMapCreate ( yMesh, WRP_Mesh, M_P_2_P, errStat2, errMsg2 )
call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
IF ( errStat >= AbortErrLev ) then
call cleanup()
RETURN
end if
CALL Transfer_Point_to_Point( yMesh, WRP_Mesh, M_P_2_P, errStat2, errMsg2, uMesh, WRP_Mesh_position ); call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
ExtBuoyancy(1:3) = WRP_Mesh%Force (:,1)
ExtBuoyancy(4:6) = WRP_Mesh%Moment(:,1)
! Write the buoyancy table headers and the external results
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A51)' ) 'Total Buoyancy loads summed about ( 0.0, 0.0, 0.0 )'
WRITE( UnSum, '(A51)' ) '---------------------------------------------------'
WRITE( UnSum, '(18x,6(2X,A20))' ) ' BuoyFxi ', ' BuoyFyi ', ' BuoyFzi ', ' BuoyMxi ', ' BuoyMyi ', ' BuoyMzi '
WRITE( UnSum, '(18x,6(2X,A20))' ) ' (N) ', ' (N) ', ' (N) ', ' (N-m) ', ' (N-m) ', ' (N-m) '
WRITE( UnSum, '(A18,6(2X,ES20.6))') 'External: ', ExtBuoyancy(1), ExtBuoyancy(2), ExtBuoyancy(3), ExtBuoyancy(4), ExtBuoyancy(5), ExtBuoyancy(6)
! Now compute internal Buoyancy
DO J = 1, yMesh%Nnodes
if (J <= numJoints) then
ptLoad = F_BF(:,J) + m%F_BF_end(:,J)
else
ptLoad = F_BF(:,J)
end if
yMesh%Force(:,J) = ptLoad(1:3)
yMesh%Moment(:,J) = ptLoad(4:6)
END DO ! DO J
IntBuoyancy = 0.0
CALL Transfer_Point_to_Point( yMesh, WRP_Mesh, M_P_2_P, errStat2, errMsg2, uMesh, WRP_Mesh_position ); call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
IntBuoyancy(1:3) = WRP_Mesh%Force(:,1)
IntBuoyancy(4:6) = WRP_Mesh%Moment(:,1)
WRITE( UnSum, '(A18,6(2X,ES20.6))') 'Internal: ', IntBuoyancy(1), IntBuoyancy(2), IntBuoyancy(3), IntBuoyancy(4), IntBuoyancy(5), IntBuoyancy(6)
IntBuoyancy = IntBuoyancy + ExtBuoyancy
WRITE( UnSum, '(A18,6(2X,ES20.6))') 'Total : ', IntBuoyancy(1), IntBuoyancy(2), IntBuoyancy(3), IntBuoyancy(4), IntBuoyancy(5), IntBuoyancy(6)
WRITE( UnSum, '(A81)') ' NOTE: External buoyancy is based on members not modelled with WAMIT'
WRITE( UnSum, '(A150)') ' Internal buoyancy is computed from all members which are marked as filled in the HydroDyn input file, regardless of PropPot flag'
WRITE( UnSum, '(A88)') ' Total buoyancy does not include WAMIT-modelled buoyancy contribution'
! ! Now compute marine growth weight at the WRP
DO J = 1, yMesh%Nnodes
if (J <= numJoints) then
yMesh%Force(:,J) = F_WMG(1:3,J) + p%F_WMG_End(:,J)
else
yMesh%Force(:,J) = F_WMG(1:3,J)
end if
yMesh%Moment(:,J) = F_WMG(4:6,J)
END DO ! DO J
MG_Wt = 0.0
CALL Transfer_Point_to_Point( yMesh, WRP_Mesh, M_P_2_P, errStat2, errMsg2, uMesh, WRP_Mesh_position ); call SetErrStat(ErrStat2, ErrMsg2, ErrStat, ErrMsg, RoutineName)
MG_Wt(1:3) = WRP_Mesh%Force(:,1)
MG_Wt(4:6) = WRP_Mesh%Moment(:,1)
!
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A36)' ) 'Weight loads about ( 0.0, 0.0, 0.0 )'
WRITE( UnSum, '(A36)' ) '------------------------------------'
WRITE( UnSum, '(18x,6(2X,A20))' ) ' MGFxi ', ' MGFyi ', ' MGFzi ', ' MGMxi ', ' MGMyi ', ' MGMzi '
WRITE( UnSum, '(18x,6(2X,A20))' ) ' (N) ', ' (N) ', ' (N) ', ' (N-m) ', ' (N-m) ', ' (N-m) '
WRITE( UnSum, '(A18,6(2X,ES20.6))') 'Marine Growth: ', MG_Wt(1), MG_Wt(2), MG_Wt(3), MG_Wt(4), MG_Wt(5), MG_Wt(6)
!
! ! Write the header for this section
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A14,I4,A44)' ) 'Nodes (first [',numJoints,'] are joints, remainder are internal nodes)'
WRITE( UnSum, '(/)' )
WRITE( UnSum, '(1X,A5,21(2X,A10))' ) ' i ', ' MbrIndx ', ' Nxi ', ' Nyi ', ' Nzi ', ' R ', ' t ', ' tMG ', ' MGDens ', ' PropPot ', 'FilledFlag', 'FilledMass', ' Cd ', ' Ca ', ' Cp ', ' Cb ', ' AxCd ', ' AxCa ', ' AxCp ', ' JAxCd ', ' JAxCa ', ' JAxCp '
WRITE( UnSum, '(1X,A5,21(2X,A10))' ) ' (-) ', ' (-) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (kg/m^3) ', ' (-) ', ' (-) ', ' (kg) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) '
! Write the node data
do I = 1,numJoints
! need to add MSL2SWL offset from this because the Positons are relative to SWL, but we should report them relative to MSL here
pos = nodes(i)%Position
pos(3) = pos(3) + p%WaveField%MSL2SWL
write( UnSum, '(1X,I5,(2X,A10),3(2X,F10.4),2(2X,A10),2(2X,ES10.3),10(2X,A10),3(2X,ES10.3))' ) i,' - ', pos, ' - ', ' - ', nodes(i)%tMG, nodes(i)%MGdensity, ' - ', ' - ', ' - ', ' - ', ' - ', ' - ', ' - ', ' - ', ' - ', ' - ', nodes(i)%JAxCd, nodes(i)%JAxCa, nodes(i)%JAxCp
end do
c = numJoints
do j= 1, numMembers
do i = 2, members(j)%NElements
c = c + 1
if (members(j)%l_fill - members(j)%dl*(i-1) > 0.0) then
fillFlag = .true.
else
fillFlag = .false.
end if
! need to add MSL2SWL offset from this because the Positons are relative to SWL, but we should report them relative to MSL here
pos = nodes(c)%Position
pos(3) = pos(3) + p%WaveField%MSL2SWL
if (members(j)%flipped) then
II=members(j)%NElements+2-I
else
II=I
endif
write( UnSum, '(1X,I5,(2X,I10),3(2X,F10.4),4(2X,ES10.3),2(6X,L6),8(2X,ES10.3),3(7x,A5))' ) c, members(j)%MemberID, pos, members(j)%R(ii), members(j)%R(ii)-members(j)%Rin(ii), members(j)%tMG(ii), members(j)%MGdensity(ii), members(j)%PropPot, fillFlag, members(j)%m_fb_u(ii)+members(j)%m_fb_l(ii), members(j)%Cd(ii), members(j)%Ca(ii), members(j)%Cp(ii), members(j)%Cb(ii), members(j)%AxCd(ii), members(j)%AxCa(ii), members(j)%AxCp(ii), ' - ', ' - ', ' - '
end do
end do
write( UnSum, '(//)' )
write( UnSum, '(A8)' ) 'Members'
write( UnSum, '(/)' )
write( UnSum, '(1X,A8,2X,A6,2X,A6,33(2X,A12))' ) 'MemberID', 'joint1','joint2',' Length ', ' NElem ', ' Volume ', ' MGVolume ', ' R1 ', ' t1 ', ' R2 ', ' t2 ', ' PropPot ', 'FilledFlag', 'FillDensity', ' FillFSLoc ', ' FillMass ', ' Cd1 ', ' Ca1 ', ' Cp1 ', ' Cb1 ', ' AxCd1 ', ' AxCa1 ', ' AxCp1 ', ' JAxCd1 ', ' JAxCa1 ', ' JAxCp1 ', ' Cd2 ', ' Ca2 ', ' Cp2 ', ' Cb2 ', ' AxCd2 ', ' AxCa2 ', ' AxCp2 ', ' JAxCd2 ', ' JAxCa2 ', ' JAxCp2 '
write( UnSum, '(1X,A8,2X,A6,2X,A6,33(2X,A12))' ) ' (-) ', ' (-) ',' (-) ',' (m) ', ' (-) ', ' (m^3) ', ' (m^3) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (-) ', ' (-) ', ' (kg/m^3) ', ' (-) ', ' (kg) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) ', ' (-) '
do i = 1,numMembers
N = members(i)%NElements
IF (members(i)%PropPot) THEN
MGvolume = 0.0
memberVol = 0.0
ELSE
memberVol = members(i)%Vouter
MGvolume = members(i)%Vouter - members(i)%Vinner
END IF
IF ( members(i)%l_fill > 0.0 ) THEN
filledFlag = .TRUE.
mass_fill = members(i)%FillDens*members(i)%Vballast
ELSE
mass_fill = 0.0
filledFlag = .FALSE.
END IF
Cd1 = members(i)%Cd(1)
Cd2 = members(i)%Cd(N+1)
Ca1 = members(i)%Ca(1)
Ca2 = members(i)%Ca(N+1)
Cp1 = members(i)%Cp(1)
Cp2 = members(i)%Cp(N+1)
AxCd1 = members(i)%AxCd(1)
AxCd2 = members(i)%AxCd(N+1)
AxCa1 = members(i)%AxCa(1)
AxCa2 = members(i)%AxCa(N+1)
AxCp1 = members(i)%AxCp(1)
AxCp2 = members(i)%AxCp(N+1)
Cb1 = members(i)%Cb(1)
Cb2 = members(i)%Cb(N+1)
JAxCd1 = nodes(members(i)%NodeIndx(1 ))%JAxCd
JAxCd2 = nodes(members(i)%NodeIndx(1+N))%JAxCd
JAxCa1 = nodes(members(i)%NodeIndx(1 ))%JAxCa
JAxCa2 = nodes(members(i)%NodeIndx(1+N))%JAxCa
JAxCp1 = nodes(members(i)%NodeIndx(1 ))%JAxCp
JAxCp2 = nodes(members(i)%NodeIndx(1+N))%JAxCp
write( UnSum, '(1X,I8,2X,I6,2X,I6,2X,ES12.5,2X,I12, 6(2X,ES12.5),2(2X,L12),23(2X,ES12.5))' ) members(i)%MemberID, &
members(i)%NodeIndx(1), members(i)%NodeIndx(N+1), members(i)%RefLength, N, &
memberVol, MGvolume, members(i)%Rmg(1), members(i)%Rmg(1)-members(i)%Rin(1), &
members(i)%Rmg(N+1), members(i)%Rmg(N+1)-members(i)%Rin(N+1), &
members(i)%PropPot, filledFlag, members(i)%FillDens, members(i)%FillFSLoc, &
mass_fill, Cd1, Ca1, Cp1, Cb1, AxCd1, AxCa1, AxCp1, JAxCd1, JAxCa1, JAxCp1, &
Cd2, Ca2, Cp2, Cb2, AxCd2, AxCa2, AxCp2, JAxCd2, JAxCa2, JAxCp2
end do ! i = 1,numMembers
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A24)' ) 'Requested Member Outputs'
WRITE( UnSum, '(/)' )
WRITE( UnSum, '(1X,A10,11(2X,A10))' ) ' Label ', ' Xi ', ' Yi ', ' Zi ', ' MemberID ', ' StartXi ', ' StartYi ', ' StartZi ', ' EndXi ', ' EndYi ', ' EndZi ', ' Loc '
WRITE( UnSum, '(1X,A10,11(2X,A10))' ) ' (-) ', ' (m) ', ' (m) ', ' (m) ', ' (-) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (m) ', ' (-) '
DO I = 1,NOutputs
tmpName = OutParam(I)%Name
IF (OutParam(I)%SignM == -1 ) tmpName = tmpName(2:10)
IF ( ( INDEX( 'mM', tmpName(1:1) ) > 0 ) .AND. (OutParam(I)%Units /= 'INVALID' ) ) THEN
!Get Member index and Node index
read (tmpName(2:2),*) mbrIndx
read (tmpName(4:4),*) nodeIndx
s = MOutLst(mbrIndx)%NodeLocs(nodeIndx)
! Find the member starting and ending node locations
! The member output is computed as a linear interpolation of the nearest two markers
mem = members(MOutLst(mbrIndx)%MemberIDIndx)
node1 = nodes(mem%NodeIndx(1))
node2 = nodes(mem%NodeIndx(mem%NElements+1))
! need to add MSL2SWL offset from this because the Positons are relative to SWL, but we should report them relative to MSL here
pos = node1%Position
pos(3) = pos(3) + p%WaveField%MSL2SWL
pos2 = node2%Position
pos2(3) = pos2(3) + p%WaveField%MSL2SWL
outLoc = pos*(1-s) + pos2*s
WRITE( UnSum, '(1X,A10,3(2x,F10.4),2x,I10,7(2x,F10.4))' ) OutParam(I)%Name, outLoc, MOutLst(mbrIndx)%MemberID, pos,pos2, s
END IF
END DO
WRITE( UnSum, '(//)' )
WRITE( UnSum, '(A24)' ) 'Requested Joint Outputs'
WRITE( UnSum, '(/)' )
WRITE( UnSum, '(1X,A10,5(2X,A10))' ) ' Label ', ' Xi ', ' Yi ', ' Zi ', 'InpJointID'
WRITE( UnSum, '(1X,A10,5(2X,A10))' ) ' (-) ', ' (m) ', ' (m) ', ' (m) ', ' (-) '
DO I = 1,NOutputs
tmpName = OutParam(I)%Name
IF (OutParam(I)%SignM == -1 ) tmpName = tmpName(2:10)
IF ( ( INDEX( 'jJ', tmpName(1:1) ) > 0 ) .AND. (OutParam(I)%Units /= 'INVALID') ) THEN
!Get Member index and Node index
read (tmpName(2:2),*) nodeIndx
m1 = JOutLst(nodeIndx)%JointIDIndx
! need to add MSL2SWL offset from this because the Positons are relative to SWL, but we should report them relative to MSL here
pos = nodes(m1)%Position
pos(3) = pos(3) + p%WaveField%MSL2SWL
WRITE( UnSum, '(1X,A10,3(2x,F10.4),2x,I10)' ) OutParam(I)%Name, pos, JOutLst(nodeIndx)%JointID
END IF
END DO
call cleanup()
contains
!...................................
subroutine cleanup()
call MeshDestroy(WRP_Mesh, ErrStat2, ErrMsg2)
call MeshDestroy(WRP_Mesh_position, ErrStat2, ErrMsg2)
call MeshMapDestroy(M_P_2_P, ErrStat2, ErrMsg2)
call Morison_DestroyNodeType(node1, ErrStat2, ErrMsg2)
call Morison_DestroyNodeType(node2, ErrStat2, ErrMsg2)
call Morison_DestroyMemberType(mem, ErrStat2, ErrMsg2)
end subroutine cleanup
END SUBROUTINE WriteSummaryFile
!----------------------------------------------------------------------------------------------------------------------------------
subroutine Morison_GenerateSimulationNodes( MSL2SWL, numJoints, inpJoints, numMembers, inpMembers, numNodes, nodes, errStat, errMsg )
! This subdivides a Morison member according to its maximum desired
! element length (MDivSize), allocating the member's arrays, and
! adding resuling new nodes to the master node array.
real(ReKi), intent (in ) :: MSL2SWL ! mean sea level To still water level offset value
integer, intent (in ) :: numJoints ! number of joints in the input file
type(Morison_JointType), intent (in ) :: inpJoints(:) ! array of input joint data structures
integer, intent (in ) :: numMembers ! number of members specified in the inputs
type(Morison_MemberInputType), intent (inout) :: inpMembers(:) ! array of input member data structures
integer, intent (inout) :: numNodes ! the total number of nodes being used for the simulation model
type(Morison_NodeType), allocatable, intent (inout) :: nodes(:) ! the array of simulation nodes
integer, intent ( out) :: errStat ! returns a non-zero value when an error occurs
character(*), intent ( out) :: errMsg ! Error message if errStat /= ErrID_None
integer :: numDiv, maxNodes
integer :: i, j
real(ReKi) :: s ! interpolation factor
real(ReKi) :: memLength ! member length
integer :: j1, j2 ! generic integer for counting
INTEGER(IntKi) :: errStat2 ! Error status of the operation (occurs after initial error)
CHARACTER(errMsgLen) :: errMsg2 ! Error message if errStat2 /= ErrID_None
! Initialize errStat
errStat = ErrID_None
errMsg = ""
! Initialize quantities
maxNodes = numJoints
! Determine maximum nodes in simulation mesh due to internal member subdivision
do i = 1,numMembers
j1 = inpMembers(I)%MJointID1Indx
j2 = inpMembers(I)%MJointID2Indx
call GetDistance(inpJoints(j1)%Position, inpJoints(j2)%Position, memLength)
if ( EqualRealNos(memLength, 0.0_ReKi) ) then
errMsg = ' Input file member with ID: '//trim(num2lstr(inpMembers(i)%MemberID))//' must have length greater than zero.'
errStat = ErrID_Fatal
return
end if
numDiv = CEILING( memLength / inpMembers(i)%MDivSize )