-
Notifications
You must be signed in to change notification settings - Fork 68
/
trees_setup_ls.f90
1862 lines (1356 loc) · 52.7 KB
/
trees_setup_ls.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
!!
!! Copyright (C) 2009-2013 Johns Hopkins University
!!
!! This file is part of lesgo.
!!
!! lesgo is free software: you can redistribute it and/or modify
!! it under the terms of the GNU General Public License as published by
!! the Free Software Foundation, either version 3 of the License, or
!! (at your option) any later version.
!!
!! lesgo is distributed in the hope that it will be useful,
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
!! GNU General Public License for more details.
!!
!! You should have received a copy of the GNU General Public License
!! along with lesgo. If not, see <http://www.gnu.org/licenses/>.
!!
module trees_setup_ls
use trees_base_ls
implicit none
save
private
public :: sdistfcn_tree_array, fill_tree_array
character (*), parameter :: mod_name = 'trees_setup'
integer :: ident = 0 !--global index used to set branch % ident
real (rp), parameter :: zero_vector(nd) = 0._rp
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!--calculates signed distance function for whole tree array
!--could use fmm to speed things up
!--when brident is present, it stores the branch % ident for pts
! inside the branches and is -1 otherwise
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine sdistfcn_tree_array (ip, mx, phi, brident)
implicit none
integer, intent (in) :: ip !--specifies chunk (or process)
! local chunk k=0 <-> ip*(mx(3)-1)
! local chunk k=nz <-> ip*(mx(3)-1)+mx(3)
integer, intent (in) :: mx(nd) ! maximum dimension in x,y,z directions
real (rp), intent (out) :: phi(mx(1), mx(2), 0:mx(3)) !--note 0 here
integer, intent (out), optional :: brident(mx(1), mx(2), mx(3))
character (*), parameter :: sub = mod_name // '.sdistfcn_tree_array'
! Make initial phi really, really (unrealistically) big
real (rp), parameter :: phi_init = huge (1._rp)
real (rp), parameter :: phi_thresh = 10._rp * epsilon (0._rp)
integer :: i
!---------------------------------------------------------------------
if (.not. grid % initialized) then
call error (sub, 'grid must be initialized')
end if
write(*,*) 'mx = ', mx
!--check size of phi array vs. size of grid
do i = 1, nd-1
if (mx(i) < grid % nx(i)) then
call error (sub, 'phi is too small in dimension', i)
end if
!--not sure how to check third dimension now with the chunks
end do
!--init phi to large values
phi = phi_init
if (present (brident)) brident = -1
do i = 1, n_tree
call mesg (sub, 'trunk % abs_dir =', tree_array(i) % trunk % abs_dir)
call mesg (sub, 'trunk % x_hat =', tree_array(i) % trunk % x_hat)
call mesg (sub, 'trunk % y_hat =', tree_array(i) % trunk % y_hat)
call mesg (sub, 'trunk % z_hat =', tree_array(i) % trunk % z_hat)
if (present (brident)) then
call sdistfcn_branch (tree_array(i) % trunk, ip, mx, phi, brident)
else
call sdistfcn_branch (tree_array(i) % trunk, ip, mx, phi)
end if
end do
!--where phi is 0 < phi < thresh, set phi to 0
!--this is just experimental--a better way must exist
!where ((0._rp < phi) .and. (phi < phi_thresh)) phi = 0._rp
where ((abs (phi) < phi_thresh)) phi = 0._rp
end subroutine sdistfcn_tree_array
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!--assumes grid is initialized, and is not larger than phi
!--puts phi on level_set_node
!--this is slow
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
recursive subroutine sdistfcn_branch (br, ip, mx, phi, brident)
implicit none
type (branch_type), intent (in) :: br
integer, intent (in) :: ip
integer, intent (in) :: mx(3)
real (rp), intent (inout) :: phi(mx(1), mx(2), 0:mx(3))
!--only 1:nx, 1:ny used
!--note 0
integer, intent (inout), optional :: brident(mx(1), mx(2), mx(3))
character (*), parameter :: sub_name = mod_name // '.sdistfcn_branch'
integer, parameter :: phi_node = 1 !--u-nodes
real (rp), parameter :: thresh = 10._rp * epsilon (0._rp)
integer :: i, j, k, ktot
real (rp) :: cosine, sine
real (rp) :: d_para, d_perp
real (rp) :: dist
real (rp) :: dist_sq
real (rp) :: l
real (rp) :: rb0, rb, rbl
real (rp) :: t, tgt
real (rp) :: x_para(nd), x_perp(nd)
real (rp) :: x(nd)
real(rp), dimension(3,3) :: lambda_skew
real(rp) :: xkeep
!---------------------------------------------------------------------
if (.not. br % resolved) return
! Length of the branch
l = br % l
t = br % taper
write(*,*) 'br%taper = ', br%taper
! Radius of the base
rb0 = 0.5_rp * (br % d)
! Radius of the top
rbl = rb0 * (1._rp - t)
tgt = rb0 * t / l
sine = rb0 * t / sqrt ((rb0 * t)**2 + l**2)
cosine = l / sqrt ((rb0 * t)**2 + l**2)
!lambda_skew=0.
!lambda_skew(1,1)=1.
!lambda_skew(2,2)=1.
!lambda_skew(3,1)=sin(45.*3.14/180);
!lambda_skew(3,3)=cos(45.*3.14/180);
do k = 0, mx(3) !--mx(3) should be (grid % nx(3) - 1) / (np) + 1
ktot = ip * (mx(3) - 1) + k
call mesg (sub_name, 'processing k =', ktot)
x(3) = pt_of_grid (ktot, 3, phi_node) - br % x0(3)
xkeep=x(3)
do j = 1, mx(2) !--mx(2) should be grid % nx(2)
x(2) = pt_of_grid (j, 2, phi_node) - br % x0(2)
do i = 1, mx(1) !--mx(1) should be grid % nx(1)
x(1) = pt_of_grid (i, 1, phi_node) - br % x0(1)
!! Perform cylinder skew operations
! x(3) = lambda_skew(3,1)*x(1) + lambda_skew(3,3)*xkeep;
!--this part depends on cross section
d_para = dot_product (x, br % abs_dir)
x_para = d_para * (br % abs_dir)
x_perp = x - x_para
d_perp = mag (x_perp)
!--set rb with taper? or explicitly code in?
rb = rb0 * (1._rp - (d_para / l) * t) !--careful when d_para < 0 or
! when d_para > l
select case (branch_cross_section)
case ('circular')
call dist_circle ()
case ('square')
call dist_square ()
case ('square+plate')
call dist_square ()
dist_sq = dist
call dist_plate ()
dist = min (dist_sq, dist)
case default
call error (sub_name, 'invalid branch_cross_section')
end select
!--need to take sign into account here
if (phi(i, j, k) >= 0._rp) then !--cannot modify interior values...
if (dist < phi(i, j, k)) then
phi(i, j, k) = dist
if (present (brident) .and. (dist <= 0._rp)) then
if (k > 0) brident(i, j, k) = br % ident
!--current convention is weird, phi from 0, brident from 1
! however, this is assumed in the main code so must change
! all or none
end if
end if
end if
end do
end do
end do
if (associated (br % sub_branch)) then
do i = 1, br % n_sub_branch
if (present (brident)) then
call sdistfcn_branch (br % sub_branch(i), ip, mx, phi, brident)
else
call sdistfcn_branch (br % sub_branch(i), ip, mx, phi)
end if
end do
end if
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
contains
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!--puts in a flat plate shaped sort of like a stealth bomber:
! *
! *****
! *********
! *************
! *****************
! ***** *****
! * *
!--this is supposed to simulate the case where the SGS branches
! basically fill in the space they occupy
!--only does something when br % gen = n_gen
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_plate ()
implicit none
real (rp) :: d1, d2, d3
real (rp) :: xi_part(nd)
real (rp) :: u(nd), u1(nd), u2(nd), u3(nd) !--temp variables
real (rp) :: xp(nd), xi(nd)
real (rp) :: xi_hat(nd), eta_hat(nd), zeta_hat(nd)
real (rp) :: hdepth, s !--dimensions of plate shape
real (rp) :: ratio
!-------------------------------------------------------------------
if (br % gen /= tree_array(br % itree) % n_gen) return
ratio = tree_array(br % itree) % ratio
hdepth = ratio * rb0 !--half-depth
s = (rb0 + (l + rb0) * ratio / (1._rp - ratio)) / sqrt (2._rp)
!--s is the half length of one of the plate sides
!--assumes x is rel. to x0 already
!--make xp the branch-local coords
xp(1) = dot_product (x, br % x_hat)
xp(2) = dot_product (x, br % y_hat)
xp(3) = dot_product (x, br % z_hat)
!--now make rel. to center of plate shape
xp = xp - l * (/ 0._rp, 0._rp, 1._rp /) !--along the branch-local z-axis
!--these are relative to branch directions (x_hat, y_hat, z_hat)
xi_hat = (/ 1._rp, 0._rp, 0._rp /) !--already normalized
eta_hat = (/ 0._rp, 1._rp, -1._rp /)
eta_hat = eta_hat / mag (eta_hat)
zeta_hat = (/ 0._rp, 1._rp, 1._rp /)
zeta_hat = zeta_hat / mag (zeta_hat)
!--rotate again to get 45 degree clockwise rotation in x_hat-z_hat plane
xi(1) = dot_product (xp, xi_hat)
xi(2) = dot_product (xp, eta_hat)
xi(3) = dot_product (xp, zeta_hat)
!--redefine xi_hat, etc. to be in their own coordinates now
xi_hat = (/ 1._rp, 0._rp, 0._rp /)
eta_hat = (/ 0._rp, 1._rp, 0._rp /)
zeta_hat = (/ 0._rp, 0._rp, 1._rp /)
!--to be subtracted from u vectors
xi_part = sign (1._rp, xi(1)) * min (abs (xi(1)), hdepth) * &
xi_hat
!--most of the inequalities here are "greedy": they include boundaries
! this ensures that points do not get missed, and is not a problem since
! the distance function is continuous across the boundaries anyway
if (xi(2) >= s .and. (0._rp <= xi(3) .and. xi(3) <= s)) then
!--case 2
u = xi - s * eta_hat - xi_part
dist = mag (u - dot_product (u, zeta_hat) * zeta_hat)
else if ((-s <= xi(2) .and. xi(2) <= 0._rp) .and. xi(3) <= -s) then
!--case 3
u = xi + s * zeta_hat - xi_part
dist = mag (u - dot_product (u, eta_hat) * eta_hat)
else if ((0._rp <= xi(2) .and. xi(2) <= s) .and. &
(-s <= xi(3) .and. xi(3) <= 0._rp)) then
!--case 4, 5
if ( abs (xi(3)) >= abs (xi(2)) ) then
!--case 4
u = xi - xi_part
dist = mag (u - dot_product (u, zeta_hat) * zeta_hat)
else
!--case 5
u = xi - xi_part
dist = mag (u - dot_product (u, eta_hat) * eta_hat)
end if
else if ((xi(2) >= 0._rp .and. xi(3) <= -s) .and. &
( abs (xi(3)) >= abs (xi(2)) )) then
!--case 6
u = xi + s * zeta_hat - xi_part
dist = mag (u) !--no proj
else if ((xi(2) >= s .and. xi(3) <= 0._rp) .and. &
( abs (xi(3)) <= abs (xi(2)) )) then
!--case 7
u = xi - s * eta_hat - xi_part
dist = mag (u) !--no proj
else if (xi(2) <= -s .and. xi(3) <= -s) then
!--case 8
u = xi + s * (eta_hat + zeta_hat) - xi_part
dist = mag (u) !--no proj
else if (xi(2) >= s .and. xi(3) >= s) then
!--case 9
u = xi - s * (eta_hat + zeta_hat) - xi_part
dist = mag (u) !--no proj
else if (xi(2) <= -s .and. abs (xi(3)) <= s) then
!--case 10
u = xi + s * (eta_hat + zeta_hat) - xi_part
dist = mag (u - dot_product (u, zeta_hat) * zeta_hat)
else if (abs (xi(2)) <= s .and. xi(3) >= s) then
!--case 11
u = xi - s * (eta_hat + zeta_hat) - xi_part
dist = mag (u - dot_product (u, eta_hat) * eta_hat)
else if (xi(2) <= -s .and. xi(3) >= s) then
!-case 12
u = xi - s * (-eta_hat + zeta_hat) - xi_part
dist = mag (u) !--no proj
else !--case 1
if (abs (xi(1)) >= hdepth) then !--outside
dist = abs (xi(1)) - hdepth
else !--inside
!--set dist after this if-block
if ((0._rp <= xi(2) .and. xi(2) <= s) .and. &
(0._rp <= xi(3) .and. xi(3) <= s)) then
!--quadrant 1
u1 = xi - xi_part
d1 = mag (u1 - dot_product (u1, eta_hat) * eta_hat)
u2 = xi - s * zeta_hat - xi_part
d2 = mag (u2 - dot_product (u2, eta_hat) * eta_hat)
u3 = xi - s * eta_hat - xi_part
d3 = mag (u3 - dot_product (u3, zeta_hat) * zeta_hat)
else if ((-s <= xi(2) .and. xi(2) <= 0._rp) .and. &
(0._rp <= xi(3) .and. xi(3) <= s)) then
!--quadrant 2
u1 = xi - xi_part
d1 = mag (u1) !--no proj
u2 = xi + s * eta_hat - xi_part
d2 = mag (u2 - dot_product (u2, zeta_hat) * zeta_hat)
u3 = xi - s * zeta_hat - xi_part
d3 = mag (u3 - dot_product (u3, eta_hat) * eta_hat)
else
!--quadrant 3
u1 = xi - xi_part
d1 = mag (u1 - dot_product (u1, zeta_hat) * zeta_hat)
u2 = xi + s * eta_hat - xi_part
d2 = mag (u2 - dot_product (u2, zeta_hat) * zeta_hat)
u3 = xi + s * zeta_hat - xi_part
d3 = mag (u3 - dot_product (u3, eta_hat) * eta_hat)
end if
dist = -min ( d1, d2, d3, hdepth - abs (xi(1)) )
end if
end if
end subroutine dist_plate
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_circle ()
implicit none
if (add_base .and. add_cap) then
call dist_circle_bc ()
else if (add_cap) then !--add_base is false
call dist_circle_c ()
else if (add_base) then !--add_cap is false
call dist_circle_b ()
else !--no bases/caps
call dist_circle_nobc ()
end if
end subroutine dist_circle
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_circle_b ()
implicit none
call error (sub_name, 'dist_circle_b not implemented')
end subroutine dist_circle_b
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_circle_c ()
implicit none
call error (sub_name, 'dist_circle_c not implemented')
end subroutine dist_circle_c
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_circle_bc ()
implicit none
!integer :: code
!-------------------------------------------------------------------
if (d_para <= 0._rp) then
if (base_shape == 'hemispherical') then
dist = mag (x) - rb0
else if (base_shape == 'rectangular') then
!--special case: allowing so-called rectangular base with
! circular cross-section
if (d_para > -rb0) then
if (-d_para > d_perp) then
dist = -(rb0 + d_para)
else
dist = d_perp - rb0
end if
else
if (d_perp < rb0) then
dist = -(rb0 + d_para)
else
dist = sqrt ((d_perp - rb0)**2 + (d_para + rb0)**2)
end if
end if
else
call error (sub_name, 'invalid base_shape')
end if
!code = 1
else if ((d_para < (d_perp - rb0) * tgt) .and. (d_perp > rb0)) then
!--know d_para > 0 here
dist = sqrt ((d_perp - rb0)**2 + d_para**2)
!code = 2
else if ((d_para < l) .and. (d_perp < rb)) then
!--know d_para > 0 here
if ((d_perp < rbl) .and. (l - d_para < (rbl - d_perp) * tgt)) then
dist = -sqrt ((l - d_para)**2 + (rbl - d_perp)**2)
!code = 3
else
dist = - (rb - d_perp) * cosine
!code = 4
end if
else if ((d_para >= l) .and. (d_perp**2 + (d_para - l)**2 <= rbl**2)) then
dist = sqrt (d_perp**2 + (d_para - l)**2) - rbl
!code = 5
else if ((d_para >= l) .and. (d_para - l >= (d_perp - rbl) * tgt)) then
dist = sqrt (d_perp**2 + (d_para - l)**2) - rbl
!code = 6
else
dist = (d_perp - rb) * cosine
!code = 7
end if
!write (*, '(1x,a,4(1x,i0),es12.5)') 'i, j, k, code, dist =', &
! i, j, k, code, dist
end subroutine dist_circle_bc
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_circle_nobc ()
implicit none
call error (sub_name, 'dist_circle_nobc not implemented yet')
end subroutine dist_circle_nobc
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_square ()
implicit none
if (add_base .and. add_cap) then
call dist_square_bc ()
else if (add_cap) then
call dist_square_c ()
else if (add_base) then
call dist_square_b ()
else
call dist_square_nobc ()
end if
end subroutine dist_square
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_square_b ()
implicit none
call error (sub_name, 'dist_square_b not implemented yet')
end subroutine dist_square_b
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_square_c ()
implicit none
call error (sub_name, 'dist_square_c not implemented yet')
end subroutine dist_square_c
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_square_bc ()
implicit none
real (rp) :: c(nd), u(nd)
real (rp) :: s1, s2
real (rp) :: x1, x2
!-------------------------------------------------------------------
x1 = dot_product (x, br % x_hat)
x2 = dot_product (x, br % y_hat)
s1 = sign (1._rp, x1)
s2 = sign (1._rp, x2)
if (d_para < -rb0) then
!--below the base
if ((abs (x1) > rb0) .and. (abs (x2) > rb0)) then
!--closest point is lower corners of base
c = rb0 * (-(br % z_hat) + s1 * (br % x_hat) + s2 * (br % y_hat))
dist = mag (x - c)
else if ((abs (x1) <= rb0) .and. (abs (x2) > rb0)) then
c = rb0 * (-(br % z_hat) + s2 * (br % y_hat))
c = x - c
c = c - (dot_product (br % x_hat, c)) * (br % x_hat)
dist = mag (c)
else if ((abs (x1) > rb0) .and. (abs (x2) <= rb0)) then
c = rb0 * (-(br % z_hat) + s1 * (br % x_hat))
c = x - c
c = c - (dot_product (br % y_hat, c)) * (br % y_hat)
dist = mag (c)
else
dist = -(d_para + rb0)
end if
else if (d_para < 0._rp) then !--know -rb0 <= d
!--level with the base
if ((abs (x1) > rb0) .and. (abs (x2) > rb0)) then
c = rb0 * (s1 * (br % x_hat) + s2 * (br % y_hat))
dist = mag (x_perp - c) !--x_perp, not x
else if ((abs (x1) <= rb0) .and. (abs (x2) > rb0)) then
dist = abs (x2) - rb0
else if ((abs (x1) > rb0) .and. (abs (x2) <= rb0)) then
dist = abs (x1) - rb0
else !--inside base
if ( abs (d_para) > max (abs (x1), abs (x2)) ) then
dist = -(rb0 - abs (d_para))
else
dist = max (abs (x1), abs (x2)) - rb0
end if
end if
else if (d_para > l + rbl) then
!--above the cap
if ((abs (x1) > rbl) .and. (abs (x2) > rbl)) then
c = (l + rbl) * (br % z_hat) + &
rbl * (s1 * (br % x_hat) + s2 * (br % y_hat))
dist = mag (x - c)
else if ((abs (x1) <= rbl) .and. (abs (x2) > rbl)) then
c = (l + rbl) * (br % z_hat) + rbl * s2 * (br % y_hat)
c = x - c
c = c - dot_product (c, br % x_hat) * (br % x_hat)
dist = mag (c)
else if ((abs (x1) > rbl) .and. (abs (x2) <= rbl)) then
c = (l + rbl) * (br % z_hat) + rbl * s1 * (br % x_hat)
c = x - c
c = c - dot_product (c, br % y_hat) * (br % y_hat)
dist = mag (c)
else
dist = d_para - (l + rbl)
end if
else if (d_para >= l) then !--know d_para <= l + rbl
if ((abs (x1) > rbl) .and. (abs (x2) > rbl)) then
c = rbl * (s1 * (br % x_hat) + s2 * (br % y_hat))
c = x_perp - c
if (d_para - l > mag (c) * tgt) then !--closer to vertical
dist = mag (c)
else !--closer to slanted line
!u = (/ -s1 * sine, -s2 * sine, cosine /)
u = (-s1 * sine) * (br % x_hat) + (-s2 * sine) * (br % y_hat) + &
cosine * (br % z_hat)
u = u / mag (u) !--unit vector along slanted corner line
c = rb0 * (s1 * (br % x_hat) + s2 * (br % y_hat))
!--reference corner slanted line passes through
c = x - c
c = c - dot_product (c, u) * u
dist = mag (c)
end if
else if ((abs (x1) <= rbl) .and. (abs (x2) > rbl)) then
if ( d_para - l > (abs (x2) - rbl) * tgt ) then !--closer to vertical
dist = abs (x2) - rbl
else !--closer to slanted plane
dist = (abs (x2) - rb) * cosine !--it is rb, not rbl!
end if
else if ((abs (x1) > rbl) .and. (abs (x2) <= rbl)) then
if ( d_para - l > (abs (x1) - rbl) * tgt ) then !--closer to vertical
dist = abs (x1) - rbl
else !--closer to slanted plane
dist = (abs (x1) - rb) * cosine !--it is rb, not rbl!
end if
else !--inside the cap
if ( d_para - l > max (abs (x1), abs (x2)) ) then
dist = -((l + rbl) - d_para)
else
dist = max (abs (x1), abs (x2)) - rbl
end if
end if
else !--d_para >= 0 & d_para <= l
if ((abs (x1) > rb) .and. (abs (x2) > rb)) then
!--check if closer to lower corner
c = rb0 * (s1 * (br % x_hat) + s2 * (br % y_hat))
if (d_para < mag (x_perp - c) * tgt) then !--corner is closer
dist = mag (x - c) !--x here, not x_perp
else !--slanted line is closer
!--this is in branch-local coordinates!
!u = (/ -s1 * sine, - s2 * sine, cosine /)
u = (-s1 * sine) * (br % x_hat) + (-s2 * sine) * (br % y_hat) + &
cosine * (br % z_hat)
u = u / mag (u)
c = rb0 * (s1 * (br % x_hat) + s2 * (br % y_hat)) !--ref. corner
c = x - c
c = c - dot_product (c, u) * u
dist = mag (c)
end if
else if ((abs (x1) <= rb) .and. (abs (x2) > rb)) then
if (d_para < (abs (x2) - rb0) * tgt) then
!--line is closer
dist = sqrt ((abs (x2) - rb0)**2 + d_para**2)
else !--slanted plane is closer
dist = (abs (x2) - rb) * cosine !--rb, not rb0
end if
else if ((abs (x1) > rb) .and. (abs (x2) <= rb)) then
if (d_para < (abs (x1) - rb0) * tgt) then
!--line is closer
dist = sqrt ((abs (x1) - rb0)**2 + d_para**2)
else !--slanted plane is closer
dist = (abs (x1) - rb) * cosine !--rb, not rb0
end if
else !--inside stem
dist = (max (abs (x1), abs (x2)) - rb) * cosine
end if
end if
end subroutine dist_square_bc
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine dist_square_nobc ()
implicit none
call error (sub_name, 'dist_square_nobc not implemented yet')
end subroutine dist_square_nobc
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
end subroutine sdistfcn_branch
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! fills the module tree array
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine fill_tree_array ()
implicit none
character (*), parameter :: sub_name = mod_name // '.fill_tree_array'
integer :: i
! for now, these parameters go here
!real (rp) :: l(n_tree) !--experimental
!real (rp) :: d(n_tree) !--experimental
!real (rp) :: x0(nd, n_tree) !--experimental
!real (rp) :: y0(n_tree)
!----------------------------------------------------------------------
!--read tree configuration file, allocates tree_array
call read_trees_conf () !--always returns n_tree > 0
!if (grid_bound) then
!
! if (grid % initialized) then ! specify in terms of grid pts
!
! !--(40, 16) for r=1/2 case (128^3)
! !--(40, 20) for r=1/2.5 case (128^3)
! !--(45, 18) for r=1/3 case (128^3)
! l = 15 * (grid % dx(3)) ! length along trunk_dir
! d = 6 * (grid % dx(1)) ! assume dx(1) = dx(2), length perp trunk_dir
!
! ! note: all three directions may be set here, but only two are
! ! actually used
! ! parallel to trunk_dir is not used here, that is set below
!
! ! just for extended case (Lx=2)
! !x0(1, 1) = ((grid % nx(1)) / 4) * (grid % dx(1))
! x0(1, 1) = ((grid % nx(1)) / 2) * (grid % dx(1))
! !x0(1, 1:2) = ((grid % nx(1)) / 4) * (grid % dx(1))
! !x0(1, 3:4) = (3 * (grid % nx(1)) / 4) * (grid % dx(1))
!
! x0(2, 1) = ((grid % nx(2)) / 2) * (grid % dx(2))
! !x0(2, (/ 1, 3 /)) = ((grid % nx(2)) / 4) * (grid % dx(2))
! !x0(2, (/ 2, 4 /)) = (3 * (grid % nx(2)) / 4) * (grid % dx(2))
!
! !x0(3, 1) = ((grid % nx(3)) / 2) * (grid % dx(3))
! x0(3, 1) = ((grid % nx(3)) / 2 - 0.5_rp) * (grid % dx(3))
! !--experiment: want to have x0 on u-node
!
! else
!
! call error (sub_name, 'for grid_bound, initialize grid first')
!
! end if
!
!else ! specify in terms of real values
!
! l = 0.2_rp
! d = 0.1_rp * l
!
! !x0(1:2) = 0.25_rp
! !x0(3:4) = 0.75_rp
!
! !y0((/ 1, 3 /)) = 0.25_rp
! !y0((/ 2, 4 /)) = 0.75_rp
!
!end if
! now make sure we start at edge of grid
! this overrides what is specified above
! old comment (still applies, but generalized to other directions now):
!! assume trees has its root at z = 0
!!tree % trunk % x0 = (/ x0, y0, 0._rp /)
!! setting x0(3) = 0. is problematic since the trunks horizontal children
!! will not lie on grid nodes.
!! experiment: shift tree up by dx(3)/2 (starts at u-node)
!! this means the trunks children should also lie on u-nodes
!call warn (sub_name, 'base of tree height may not be set correctly')
!i = maxloc (abs (trunk_dir), dim=1)
!
!if (grid_bound .and. grid % initialized) then
!
! x0(i, :) = grid % x_min(i, tree_node)
!
!else
!
! call warn (sub_name, 'grid not initialized, ' // &
! 'setting base of tree at height 0')
! x0(i, :) = 0._rp
!
!end if
do i = 1, n_tree
call init_tree (i, tree_array(i))
end do
end subroutine fill_tree_array
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!--provides a valid n_tree > 0
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine read_trees_conf ()
use messages
use string_util, only : eat_whitespace
implicit none
character (*), parameter :: sub = 'read_trees_conf'
character (*), parameter :: ftrees_conf = 'trees.conf'
character (*), parameter :: comment = '#'
character (*), parameter :: ldelim = '(' !--no whitespace allowed
character (*), parameter :: rdelim = ')' !--no whitespace allowed
character (*), parameter :: begin_tree_block = '{'
character (*), parameter :: end_tree_block = '}'
character (*), parameter :: esyntax = 'syntax error at line'
integer, parameter :: lun = 1
integer, parameter :: BUFF_LEN = 256
character (BUFF_LEN) :: buff
integer :: eqpos
integer :: i
integer :: i_tree
integer :: ios
integer :: line
integer :: n_sub_branch
integer :: state
logical :: exst
logical :: have_n_tree, have_n_sub_branch, have_rel_dir, have_root_height
!---------------------------------------------------------------------
inquire (file=ftrees_conf, exist=exst)
if (exst) then
open (lun, file=ftrees_conf, action='read')
else
call error (sub, 'file ' // ftrees_conf // ' does not exist')
end if
line = 0
have_n_tree = .false.
have_n_sub_branch = .false.
state = 0
i_tree = 0
!--format description
! * 'comment' is comment token, must be first nonblank character
! * blank lines are allowed
! * <variable name> = <variable value> set the variables: e.g.
! n_tree = 4
! * trees are defined using tree blocks of the form
! tree = {
! l = ...
! d = ...
! x0 = ...
! ...
! }
! the {,} delimiters are reqd to be in the lines as shown
do
read (lun, '(a)', iostat=ios) buff
if (ios /= 0) exit
line = line + 1
!--remove leading/intermediate whitespace
call eat_whitespace (buff)
if (verify (buff, ' ') == 0) cycle !--drop blank lines