SABER
tools_stripack.F90
Go to the documentation of this file.
1 # 1 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/tools_stripack.fypp"
2 # 1 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/../instrumentation.fypp" 1
3 # 1 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/../subr_list.fypp" 1
4 !----------------------------------------------------------------------
5 ! Header: subr_list
6 !> Subroutines/functions list
7 ! Author: Benjamin Menetrier
8 ! Licensing: this code is distributed under the CeCILL-C license
9 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
10 !----------------------------------------------------------------------
11 
12 # 926 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/../subr_list.fypp"
13 # 2 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/../instrumentation.fypp" 2
14 !----------------------------------------------------------------------
15 ! Header: instrumentation
16 !> Instrumentation functions
17 ! Author: Benjamin Menetrier
18 ! Licensing: this code is distributed under the CeCILL-C license
19 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
20 !----------------------------------------------------------------------
21 
22 # 112 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/../instrumentation.fypp"
23 # 2 "/Users/miesch/JEDI/code/working_copy/internal/mpas-bundle/saber/src/saber/external/tools_stripack.fypp" 2
24 !----------------------------------------------------------------------
25 ! Module: tools_stripack
26 !> STRIPACK routines
27 ! Source: https://people.sc.fsu.edu/~jburkardt/f_src/stripack/stripack.html
28 ! Author: Robert Renka
29 ! Original licensing: none
30 ! Modified by Benjamin Menetrier for BUMP
31 ! Licensing: this code is distributed under the CeCILL-C license
32 ! Copyright 2015-... UCAR, CERFACS, METEO-FRANCE and IRIT
33 !----------------------------------------------------------------------
35 
36 use tools_const, only: zero,one,four,hundred,pi
37 use tools_kinds, only: kind_real
38 use tools_repro, only: eq,inf,sup,supeq,small
39 use type_mpl, only: mpl_type
40 
41 
42 implicit none
43 
44 interface addnod
45  module procedure stripack_addnod
46 end interface
47 interface area
48  module procedure stripack_area
49 end interface
50 interface bdyadd
51  module procedure stripack_bdyadd
52 end interface
53 interface bnodes
54  module procedure stripack_bnodes
55 end interface
56 interface covsph
57  module procedure stripack_covsph
58 end interface
59 interface det
60  module procedure stripack_det
61 end interface
62 interface insert
63  module procedure stripack_insert
64 end interface
65 interface intadd
66  module procedure stripack_intadd
67 end interface
68 interface jrand
69  module procedure stripack_jrand
70 end interface
71 interface left
72  module procedure stripack_left
73 end interface
74 interface lstptr
75  module procedure stripack_lstptr
76 end interface
77 interface swap
78  module procedure stripack_swap
79 end interface
80 interface swptst
81  module procedure stripack_swptst
82 end interface
83 interface trfind
84  module procedure stripack_trfind
85 end interface
86 interface trmesh
87  module procedure stripack_trmesh
88 end interface
89 
90 private
91 public :: area,bnodes,trfind,trmesh
92 
93 contains
94 
95 !----------------------------------------------------------------------
96 ! Subroutine: stripack_addnod
97 !> Add a node to a triangulation
98 !----------------------------------------------------------------------
99 subroutine stripack_addnod(mpl,nst,k,x,y,z,list,lptr,lend,lnew)
100 
101 implicit none
102 
103 ! Passed variables
104 type(mpl_type),intent(inout) :: mpl !< MPI data
105 integer,intent(in) :: nst !< Index of a node at which TRFIND, begins its search.Search time depends on the proximity of this node to K. If NST<1, the search is begun at node K-1.
106 integer,intent(in) :: k !< Nodal index (index for X, Y, Z, and LEND) of the new node to be added. 4<=K.
107 real(kind_real),intent(in) :: x(k) !< X-coordinates of the nodes
108 real(kind_real),intent(in) :: y(k) !< Y-coordinates of the nodes
109 real(kind_real),intent(in) :: z(k) !< Z-coordinates of the nodes
110 integer,intent(inout) :: list(*) !< On input, the data structure associated with the triangulation of nodes 1 to K-1. On output, the data has been updated to include node K. The array lengths are assumed to be large enough to add node K. Refer to TRMESH.
111 integer,intent(inout) :: lptr(*) !< cf. list
112 integer,intent(inout) :: lend(k) !< cf. list
113 integer,intent(inout) :: lnew !< cf. list
114 
115 ! Local variables
116 integer :: i1,i2,i3,in1,io1,io2,ist,kk,km1,l,lp,lpf,lpo1,lpo1s
117 real(kind_real) :: b1,b2,b3,p(3)
118 logical :: output
119 
120 ! Set name
121 
122 
123 ! Probe in
124 
125 
126 ! Initialization
127 kk = k
128 if (kk<4) call mpl%abort('stripack_addnod','K<4')
129 km1 = kk-1
130 ist = nst
131 if (ist<1) ist = km1
132 p(1) = x(kk)
133 p(2) = y(kk)
134 p(3) = z(kk)
135 
136 ! Find a triangle (I1,I2,I3) containing K or the rightmost (I1) and leftmost (I2) visible boundary nodes as viewed from node K.
137 call trfind(ist,p,km1,x,y,z,list,lptr,lend,b1,b2,b3,i1,i2,i3)
138 
139 ! Test for colinear or duplicate nodes.
140 if (i1==0) call mpl%abort('stripack_addnod','colinear points')
141 if (i3/=0) then
142  l = i1
143  if (eq(p(1),x(l)).and.eq(p(2),y(l)).and.eq(p(3),z(l))) call mpl%abort('stripack_addnod','colinear or duplicate nodes')
144  l = i2
145  if (eq(p(1),x(l)).and.eq(p(2),y(l)).and.eq(p(3),z(l))) call mpl%abort('stripack_addnod','colinear or duplicate nodes')
146  l = i3
147  if (eq(p(1),x(l)).and.eq(p(2),y(l)).and.eq(p(3),z(l))) call mpl%abort('stripack_addnod','colinear or duplicate nodes')
148  call intadd(kk,i1,i2,i3,list,lptr,lend,lnew)
149 else
150  if (i1/=i2) then
151  call bdyadd(kk,i1,i2,list,lptr,lend,lnew)
152  else
153  call covsph(kk,i1,list,lptr,lend,lnew)
154  end if
155 end if
156 
157 ! Initialize variables for optimization of the triangulation.
158 lp = lend(kk)
159 lpf = lptr(lp)
160 io2 = list(lpf)
161 lpo1 = lptr(lpf)
162 io1 = abs(list(lpo1))
163 
164 ! Begin loop: find the node opposite K.
165 do
166  call lstptr(lend(io1),io2,list,lptr,lp)
167  if (0<=list(lp)) then
168  lp = lptr(lp)
169  in1 = abs(list(lp))
170 
171  ! Swap test: if a swap occurs,two new arcs are opposite K and must be tested.
172  lpo1s = lpo1
173  call swptst(in1,kk,io1,io2,x,y,z,output)
174  if (.not.output) then
175  if ((lpo1==lpf).or.(list(lpo1)<0)) exit
176  io2 = io1
177  lpo1 = lptr(lpo1)
178  io1 = abs(list(lpo1))
179  cycle
180  end if
181  call swap(in1,kk,io1,io2,list,lptr,lend,lpo1)
182 
183  ! A swap is not possible because KK and IN1 are already adjacent. This error in SWPTST only occurs in the neutral case and when there are nearly duplicate nodes.
184  if (lpo1/=0) then
185  io1 = in1
186  cycle
187  end if
188  lpo1 = lpo1s
189  end if
190 
191  ! No swap occurred. Test for termination and reset IO2 and IO1.
192  if ((lpo1==lpf).or.(list(lpo1)<0)) exit
193  io2 = io1
194  lpo1 = lptr(lpo1)
195  io1 = abs(list(lpo1))
196 end do
197 
198 ! Probe out
199 
200 
201 end subroutine stripack_addnod
202 
203 !----------------------------------------------------------------------
204 ! Subroutine: stripack_area
205 !> Compute the area of a spherical triangle
206 !----------------------------------------------------------------------
207 subroutine stripack_area(v1,v2,v3,area)
208 
209 implicit none
210 
211 ! Passed variables
212 real(kind_real),intent(in) :: v1(3) !< First point cartesian coordinates
213 real(kind_real),intent(in) :: v2(3) !< Second point cartesian coordinates
214 real(kind_real),intent(in) :: v3(3) !< Third point cartesian coordinates
215 real(kind_real),intent(out) :: area !< Area on the unit sphere
216 
217 ! Local variables
218 real(kind_real) :: a1,a2,a3,ca1,ca2,ca3,s12,s23,s31
219 real(kind_real) :: dv1(3),dv2(3),dv3(3),u12(3),u23(3),u31(3)
220 
221 ! Set name
222 
223 
224 ! Probe in
225 
226 
227 ! Initialization
228 dv1 = v1
229 dv2 = v2
230 dv3 = v3
231 
232 ! Compute cross products Uij = Vi X Vj
233 u12(1) = dv1(2)*dv2(3)-dv1(3)*dv2(2)
234 u12(2) = dv1(3)*dv2(1)-dv1(1)*dv2(3)
235 u12(3) = dv1(1)*dv2(2)-dv1(2)*dv2(1)
236 u23(1) = dv2(2)*dv3(3)-dv2(3)*dv3(2)
237 u23(2) = dv2(3)*dv3(1)-dv2(1)*dv3(3)
238 u23(3) = dv2(1)*dv3(2)-dv2(2)*dv3(1)
239 u31(1) = dv3(2)*dv1(3)-dv3(3)*dv1(2)
240 u31(2) = dv3(3)*dv1(1)-dv3(1)*dv1(3)
241 u31(3) = dv3(1)*dv1(2)-dv3(2)*dv1(1)
242 
243 ! Normalize Uij to unit vectors
244 s12 = dot_product(u12,u12)
245 s23 = dot_product(u23,u23)
246 s31 = dot_product(u31,u31)
247 
248 ! Test for a degenerate triangle associated with collinear vertices
249 if ((.not.(abs(s12)>zero)).or.(.not.(abs(s23)>zero)).or.(.not.(abs(s31)>zero))) then
250  area = zero
251 else
252  s12 = sqrt(s12)
253  s23 = sqrt(s23)
254  s31 = sqrt(s31)
255  u12 = u12/s12
256  u23 = u23/s23
257  u31 = u31/s31
258 
259  ! Compute interior angles Ai as the dihedral angles between planes
260  ca1 = -dot_product(u12,u31)
261  ca2 = -dot_product(u23,u12)
262  ca3 = -dot_product(u31,u23)
263  ca1 = max(ca1,-one)
264  ca1 = min(ca1,+one)
265  ca2 = max(ca2,-one)
266  ca2 = min(ca2,+one)
267  ca3 = max(ca3,-one)
268  ca3 = min(ca3,+one)
269  a1 = acos(ca1)
270  a2 = acos(ca2)
271  a3 = acos(ca3)
272 
273  ! Compute area = a1 + a2 + a3 - pi
274  area = a1+a2+a3-pi
275  if (area<zero) area = zero
276 end if
277 
278 ! Probe out
279 
280 
281 end subroutine stripack_area
282 
283 
284 !----------------------------------------------------------------------
285 ! Subroutine: stripack_bdyadd
286 !> Add a boundary node to a triangulation
287 !----------------------------------------------------------------------
288 subroutine stripack_bdyadd(kk,i1,i2,list,lptr,lend,lnew)
289 
290 implicit none
291 
292 ! Passed variable
293 integer,intent(in) :: kk !< Index of a node to be connected to the sequence of all visible boundary nodes. 1<=KK and KK must not be equal to I1 or I2.
294 integer,intent(in) :: i1 !< First (rightmost as viewed from KK) boundary node in the triangulation that is visible from node KK (the line segment KK-I1 intersects no arcs.
295 integer,intent(in) :: i2 !< Last (leftmost) boundary node that is visible from node KK. I1 and I2 may be determined by TRFIND.
296 integer,intent(inout) :: list(*) !< Triangulation data structure created by TRMESH. Nodes I1 and I2 must be included in the triangulation. On output, the data structure is updated with the addition of node KK. Node KK is connected to I1, I2, and all boundary nodes in between.
297 integer,intent(inout) :: lptr(*) !< cf. list
298 integer,intent(inout) :: lend(*) !< cf. list
299 integer,intent(inout) :: lnew !< cf. list
300 
301 ! Local variables
302 integer :: k,lp,lsav,n1,n2,next,nsav
303 
304 ! Set name
305 
306 
307 ! Probe in
308 
309 
310 ! Initialization
311 k = kk
312 n1 = i1
313 n2 = i2
314 
315 ! Add K as the last neighbor of N1.
316 lp = lend(n1)
317 lsav = lptr(lp)
318 lptr(lp) = lnew
319 list(lnew) = -k
320 lptr(lnew) = lsav
321 lend(n1) = lnew
322 lnew = lnew+1
323 next = -list(lp)
324 list(lp) = next
325 nsav = next
326 
327 ! Loop on the remaining boundary nodes between N1 and N2, adding K as the first neighbor.
328 do
329  lp = lend(next)
330  call insert(k,lp,list,lptr,lnew)
331  if (next==n2) exit
332  next = -list(lp)
333  list(lp) = next
334 end do
335 
336 ! Add the boundary nodes between N1 and N2 as neighbors of node K.
337 lsav = lnew
338 list(lnew) = n1
339 lptr(lnew) = lnew+1
340 lnew = lnew+1
341 next = nsav
342 do
343  if (next==n2) exit
344  list(lnew) = next
345  lptr(lnew) = lnew+1
346  lnew = lnew+1
347  lp = lend(next)
348  next = list(lp)
349 end do
350 list(lnew) = -n2
351 lptr(lnew) = lsav
352 lend(k) = lnew
353 lnew = lnew+1
354 
355 ! Probe out
356 
357 
358 end subroutine stripack_bdyadd
359 
360 !----------------------------------------------------------------------
361 ! Subroutine: stripack_bnodes
362 !> Return the boundary nodes of a triangulation
363 !----------------------------------------------------------------------
364 subroutine stripack_bnodes(n,list,lptr,lend,nodes,nb)
365 
366 implicit none
367 
368 ! Passed variables
369 integer,intent(in) :: n !< Number of nodes in the triangulation. 3 <= N.
370 integer,intent(in) :: list(6*(n-2)) !< The data structure defining the triangulation, created by TRMESH.
371 integer,intent(in) :: lptr(6*(n-2)) !< cf. list
372 integer,intent(in) :: lend(n) !< cf. list
373 integer,intent(out) :: nodes(n) !< Ordered sequence of NB boundary node indexes in the range 1 to N.
374 integer,intent(out) :: nb !< Number of boundary nodes.
375 
376 ! Local variables
377 integer :: i,k,lp,n0,na,nn,nst,nt
378 
379 ! Set name
380 
381 
382 ! Probe in
383 
384 
385 ! Initialization
386 nn = n
387 
388 ! Search for a boundary node
389 nst = 0
390 do i = 1, nn
391  lp = lend(i)
392  if (list(lp)<0) then
393  nst = i
394  exit
395  end if
396 end do
397 
398 ! The triangulation contains no boundary nodes
399 if (nst==0) then
400  nb = 0
401  na = 3*(nn-2)
402  nt = 2*(nn-2)
403 
404  return
405 end if
406 
407 ! NST is the first boundary node encountered
408 ! Initialize for traversal of the boundary
409 nodes(1) = nst
410 k = 1
411 n0 = nst
412 
413 ! Traverse the boundary in counterclockwise order
414 do
415  lp = lend(n0)
416  lp = lptr(lp)
417  n0 = list(lp)
418  if (n0==nst) then
419  exit
420  end if
421  k = k+1
422  nodes(k) = n0
423 end do
424 
425 ! Store the counts
426 nb = k
427 nt = 2*n-nb-2
428 na = nt+n-1
429 
430 ! Probe out
431 
432 
433 end subroutine stripack_bnodes
434 
435 !----------------------------------------------------------------------
436 ! Subroutine: stripack_covsph
437 !> Connect an exterior node to boundary nodes, covering the sphere
438 !----------------------------------------------------------------------
439 subroutine stripack_covsph(kk,n0,list,lptr,lend,lnew)
440 
441 implicit none
442 
443 ! Passed variables
444 integer,intent(in) :: kk !< Index of the node to be connected to the set of all boundary nodes. 4<=KK.
445 integer,intent(in) :: n0 !< Index of a boundary node (in the range 1 to KK-1). N0 may be determined by TRFIND.
446 integer,intent(inout) :: list(*) !< Triangulation data structure created by TRMESH. Node N0 must be included in the triangulation. On output, updated with the addition of node KK as the last entry. The updated triangulation contains no boundary nodes.
447 integer,intent(inout) :: lptr(*) !< cf. list
448 integer,intent(inout) :: lend(*) !< cf. list
449 integer,intent(inout) :: lnew !< cf. list
450 
451 ! Local variables
452 integer :: k,lp,lsav,next,nst
453 
454 ! Set name
455 
456 
457 ! Probe in
458 
459 
460 ! Initialization
461 k = kk
462 nst = n0
463 
464 ! Traverse the boundary in clockwise order, inserting K as the first neighbor of each boundary node, and converting the boundary node to an interior node.
465 next = nst
466 do
467  lp = lend(next)
468  call insert(k,lp,list,lptr,lnew)
469  next = -list(lp)
470  list(lp) = next
471  if (next==nst) exit
472 end do
473 
474 ! Traverse the boundary again, adding each node to K's adjacency list.
475 lsav = lnew
476 do
477  lp = lend(next)
478  list(lnew) = next
479  lptr(lnew) = lnew+1
480  lnew = lnew+1
481  next = list(lp)
482  if (next==nst) exit
483 end do
484 lptr(lnew-1) = lsav
485 lend(k) = lnew-1
486 
487 ! Probe out
488 
489 
490 end subroutine stripack_covsph
491 
492 !----------------------------------------------------------------------
493 ! Subroutine: stripack_det
494 !> Compute 3D determinant
495 !----------------------------------------------------------------------
496 subroutine stripack_det(x1,y1,z1,x2,y2,z2,x0,y0,z0,output)
497 
498 implicit none
499 
500 ! Passed variables
501 real(kind_real),intent(in) :: x1 !< X-coordinate, term 1
502 real(kind_real),intent(in) :: y1 !< Y-coordinate, term 1
503 real(kind_real),intent(in) :: z1 !< Z-coordinate, term 1
504 real(kind_real),intent(in) :: x2 !< X-coordinate, term 2
505 real(kind_real),intent(in) :: y2 !< Y-coordinate, term 2
506 real(kind_real),intent(in) :: z2 !< Z-coordinate, term 2
507 real(kind_real),intent(in) :: x0 !< X-coordinate, term 0
508 real(kind_real),intent(in) :: y0 !< Y-coordinate, term 0
509 real(kind_real),intent(in) :: z0 !< Z-coordinate, term 0
510 real(kind_real),intent(out) :: output !< Determinant
511 
512 ! Local variables
513 real(kind_real) :: t1,t2,t3
514 
515 ! Set name
516 
517 
518 ! Probe in
519 
520 
521 ! Compute determinant terms
522 t1 = x0*(y1*z2-y2*z1)
523 t2 = y0*(x1*z2-x2*z1)
524 t3 = z0*(x1*y2-x2*y1)
525 
526 ! Determinant
527 output = t1-t2+t3
528 
529 ! Indistinguishability threshold for cross-platform reproducibility
530 if (small(output,t1).or.small(output,t2).or.small(output,t3)) output = zero
531 
532 ! Probe out
533 
534 
535 end subroutine stripack_det
536 
537 !----------------------------------------------------------------------
538 ! Subroutine: stripack_insert
539 !> Insert K as a neighbor of N1
540 !----------------------------------------------------------------------
541 subroutine stripack_insert(k,lp,list,lptr,lnew)
542 
543 implicit none
544 
545 ! Passed variables
546 integer,intent(in) :: k !< Index of the node to be inserted.
547 integer,intent(in) :: lp !< LIST pointer of N2 as a neighbor of N1.
548 integer,intent(inout) :: list(*) !< Triangulation data structure created by TRMESH. On output, updated with the addition of node K.
549 integer,intent(inout) :: lptr(*) !< cf. list
550 integer,intent(inout) :: lnew !< cf. list
551 
552 ! Local variables
553 integer :: lsav
554 
555 ! Set name
556 
557 
558 ! Probe in
559 
560 
561 lsav = lptr(lp)
562 lptr(lp) = lnew
563 list(lnew) = k
564 lptr(lnew) = lsav
565 lnew = lnew+1
566 
567 ! Probe out
568 
569 
570 end subroutine stripack_insert
571 
572 !----------------------------------------------------------------------
573 ! Subroutine: stripack_intadd
574 !> Add an interior node to a triangulation
575 !----------------------------------------------------------------------
576 subroutine stripack_intadd(kk,i1,i2,i3,list,lptr,lend,lnew)
577 
578 implicit none
579 
580 ! Passed variables
581 integer,intent(in) :: kk !< Index of the node to be inserted. 1<=KK and KK must not be equal to I1, I2, or I3.
582 integer,intent(in) :: i1 !< First index of the counterclockwise-ordered sequence of vertices of a triangle which contains node KK.
583 integer,intent(in) :: i2 !< Second index.
584 integer,intent(in) :: i3 !< Third index.
585 integer,intent(inout) :: list(*) !< Triangulation data structure created by TRMESH. Triangle (I1,I2,I3) must be included in the triangulation. On output, updated with the addition of node KK. KK will be connected to nodes I1, I2, and I3.
586 integer,intent(inout) :: lptr(*) !< cf. list
587 integer,intent(inout) :: lend(*) !< cf. list
588 integer,intent(inout) :: lnew !< cf. list
589 
590 ! Local variables
591 integer :: k,lp,n1,n2,n3
592 
593 ! Set name
594 
595 
596 ! Probe in
597 
598 
599 ! Initialization.
600 k = kk
601 n1 = i1
602 n2 = i2
603 n3 = i3
604 
605 ! Add K as a neighbor of I1, I2, and I3.
606 call lstptr(lend(n1),n2,list,lptr,lp)
607 call insert(k,lp,list,lptr,lnew)
608 call lstptr(lend(n2),n3,list,lptr,lp)
609 call insert(k,lp,list,lptr,lnew)
610 call lstptr(lend(n3),n1,list,lptr,lp)
611 call insert(k,lp,list,lptr,lnew)
612 
613 ! Add I1, I2, and I3 as neighbors of K.
614 list(lnew) = n1
615 list(lnew+1) = n2
616 list(lnew+2) = n3
617 lptr(lnew) = lnew+1
618 lptr(lnew+1) = lnew+2
619 lptr(lnew+2) = lnew
620 lend(k) = lnew+2
621 lnew = lnew+3
622 
623 ! Probe out
624 
625 
626 end subroutine stripack_intadd
627 
628 !----------------------------------------------------------------------
629 ! Subroutine: stripack_jrand
630 !> Return a random integer between 1 and N
631 !----------------------------------------------------------------------
632 subroutine stripack_jrand(n,ix,iy,iz,output)
633 
634 implicit none
635 
636 ! Passed variables
637 integer,intent(in) :: n !< Maximum value to be returned.
638 integer,intent(inout) :: ix !< First seed initialized to values in the range 1 to 30,000 before the first call to JRAND, and not altered between subsequent calls (unless a sequence of random numbers is to be repeated by reinitializing the seeds).
639 integer,intent(inout) :: iy !< Second seed.
640 integer,intent(inout) :: iz !< Third seed.
641 integer,intent(out) :: output !< A random integer in the range 1 to N.
642 
643 ! Local variables
644 real(kind_real) :: u,x
645 
646 ! Set name
647 
648 
649 ! Probe in
650 
651 
652 ! Initialization
653 ix = mod(171*ix,30269)
654 iy = mod(172*iy,30307)
655 iz = mod(170*iz,30323)
656 
657 ! Get x
658 x = (real(ix,kind_real)/30269.0_kind_real) &
659  & + (real(iy,kind_real)/30307.0_kind_real) &
660  & + (real(iz,kind_real)/30323.0_kind_real)
661 
662 ! Get u
663 u = x-int(x)
664 
665 ! Random integer
666 output = int(real(n,kind_real)*u)+1
667 
668 ! Probe out
669 
670 
671 end subroutine stripack_jrand
672 
673 !----------------------------------------------------------------------
674 ! Subroutine: stripack_left
675 !> Determine whether a node is to the left of a plane through the origin
676 !----------------------------------------------------------------------
677 subroutine stripack_left(x1,y1,z1,x2,y2,z2,x0,y0,z0,output)
678 
679 implicit none
680 
681 ! Passed variables
682 real(kind_real),intent(in) :: x1 !< X-coordinate, term 1
683 real(kind_real),intent(in) :: y1 !< Y-coordinate, term 1
684 real(kind_real),intent(in) :: z1 !< Z-coordinate, term 1
685 real(kind_real),intent(in) :: x2 !< X-coordinate, term 2
686 real(kind_real),intent(in) :: y2 !< Y-coordinate, term 2
687 real(kind_real),intent(in) :: z2 !< Z-coordinate, term 2
688 real(kind_real),intent(in) :: x0 !< X-coordinate, term 0
689 real(kind_real),intent(in) :: y0 !< Y-coordinate, term 0
690 real(kind_real),intent(in) :: z0 !< Z-coordinate, term 0
691 logical,intent(out) :: output !< TRUE if and only if N0 is in the closed left hemisphere.
692 
693 ! Local variables
694 real(kind_real) :: zz
695 
696 ! Set name
697 
698 
699 ! Probe in
700 
701 
702 ! LEFT = TRUE iff <N0,N1 X N2> = det(N0,N1,N2) >= 0.
703 call det(x1,y1,z1,x2,y2,z2,x0,y0,z0,zz)
704 output = sup(zz,zero)
705 
706 ! Probe out
707 
708 
709 end subroutine stripack_left
710 
711 !----------------------------------------------------------------------
712 ! Subroutine: stripack_lstptr
713 !> Return the index of NB in the adjacency list
714 !----------------------------------------------------------------------
715 subroutine stripack_lstptr(lpl,nb,list,lptr,output)
716 
717 implicit none
718 
719 ! Passed variables
720 integer,intent(in) :: lpl !< Equal to LEND(N0).
721 integer,intent(in) :: nb !< Index of the node whose pointer is to be returned. NB must be connected to N0.
722 integer,intent(in) :: list(*) !< Triangulation data structure created by TRMESH. Triangle (I1,I2,I3) must be included in the triangulation. On output, updated with the addition of node KK. KK will be connected to nodes I1, I2, and I3.
723 integer,intent(in) :: lptr(*) !< cf. list
724 integer,intent(out) :: output !< Pointer such that LIST(output) = NB or LIST(output) = -NB, unless NB is not a neighbor of N0, in which case output = LPL.
725 
726 ! Local variables
727 integer :: lp,nd
728 
729 ! Set name
730 
731 
732 ! Probe in
733 
734 
735 ! Initialization
736 lp = lptr(lpl)
737 
738 do
739  nd = list(lp)
740  if (nd==nb) exit
741  lp = lptr(lp)
742  if (lp==lpl) exit
743 end do
744 output = lp
745 
746 ! Probe out
747 
748 
749 end subroutine stripack_lstptr
750 
751 !----------------------------------------------------------------------
752 ! Subroutine: stripack_swap
753 !> Replace the diagonal arc of a quadrilateral with the other diagonal
754 !----------------------------------------------------------------------
755 subroutine stripack_swap(in1,in2,io1,io2,list,lptr,lend,lp21)
756 
757 implicit none
758 
759 ! Passed variables
760 integer,intent(in) :: in1 !< First nodal index of the vertices of the quadrilateral. IO1-IO2 is replaced by IN1-IN2. (IO1,IO2,IN1) and (IO2,IO1,IN2) must be triangles on input.
761 integer,intent(in) :: in2 !< Second nodal index.
762 integer,intent(in) :: io1 !< Third nodal index.
763 integer,intent(in) :: io2 !< Fourth nodal index.
764 integer,intent(inout) :: list(*) !< Triangulation data structure created by TRMESH. On output, updated with the swap; triangles (IO1,IO2,IN1) an (IO2,IO1,IN2) are replaced by (IN1,IN2,IO2) and (IN2,IN1,IO1) unless LP21 = 0.
765 integer,intent(inout) :: lptr(*) !< cf. list
766 integer,intent(inout) :: lend(*) !< cf. list
767 integer,intent(out) :: lp21 !< Index of IN1 as a neighbor of IN2 after the swap is performed unless IN1 and IN2 are adjacent on input, in which case LP21 = 0.
768 
769 ! Local variables
770 integer :: lp,lph,lpsav
771 
772 ! Set name
773 
774 
775 ! Probe in
776 
777 
778 ! Test for IN1 and IN2 adjacent.
779 call lstptr(lend(in1),in2,list,lptr,lp)
780 if (abs(list(lp))==in2) then
781  lp21 = 0
782 
783  return
784 end if
785 
786 ! Delete IO2 as a neighbor of IO1.
787 call lstptr(lend(io1),in2,list,lptr,lp)
788 lph = lptr(lp)
789 lptr(lp) = lptr(lph)
790 
791 ! If IO2 is the last neighbor of IO1, make IN2 the last neighbor.
792 if (lend(io1)==lph) lend(io1) = lp
793 
794 ! Insert IN2 as a neighbor of IN1 following IO1 using the hole created above.
795 call lstptr(lend(in1),io1,list,lptr,lp)
796 lpsav = lptr(lp)
797 lptr(lp) = lph
798 list(lph) = in2
799 lptr(lph) = lpsav
800 
801 ! Delete IO1 as a neighbor of IO2.
802 call lstptr(lend(io2),in1,list,lptr,lp)
803 lph = lptr(lp)
804 lptr(lp) = lptr(lph)
805 
806 ! If IO1 is the last neighbor of IO2, make IN1 the last neighbor.
807 if (lend(io2)==lph) lend(io2) = lp
808 
809 ! Insert IN1 as a neighbor of IN2 following IO2.
810 call lstptr(lend(in2),io2,list,lptr,lp)
811 lpsav = lptr(lp)
812 lptr(lp) = lph
813 list(lph) = in1
814 lptr(lph) = lpsav
815 lp21 = lph
816 
817 ! Probe out
818 
819 
820 end subroutine stripack_swap
821 
822 !----------------------------------------------------------------------
823 ! Subroutine: stripack_swptst
824 !> Decide whether to replace a diagonal arc by the other
825 !----------------------------------------------------------------------
826 subroutine stripack_swptst(n1,n2,n3,n4,x,y,z,output)
827 
828 implicit none
829 
830 integer,intent(in) :: n1 !< First index of the four nodes defining the quadrilateral with N1 adjacent to N2, and (N1,N2,N3) in counterclockwise order. The arc connecting N1 to N2 should be replaced by an arc connecting N3 to N4 if SWPTST = TRUE. Refer to subroutine SWAP.
831 integer,intent(in) :: n2 !< Second index.
832 integer,intent(in) :: n3 !< Third index.
833 integer,intent(in) :: n4 !< Fourth index.
834 real(kind_real),intent(in) :: x(*) !< X-coordinate of the nodes.
835 real(kind_real),intent(in) :: y(*) !< Y-coordinate of the nodes.
836 real(kind_real),intent(in) :: z(*) !< Z-coordinate of the nodes.
837 logical,intent(out) :: output !< TRUE if and only if the arc connecting N1 and N2 should be swapped for an arc connecting N3 and N4.
838 
839 ! Local variables
840 real(kind_real) :: dx1,dx2,dx3,dy1,dy2,dy3,dz1,dz2,dz3,x4,y4,z4,zz
841 
842 ! Set name
843 
844 
845 ! Probe in
846 
847 
848 ! Initialization
849 x4 = x(n4)
850 y4 = y(n4)
851 z4 = z(n4)
852 dx1 = x(n1)-x4
853 dx2 = x(n2)-x4
854 dx3 = x(n3)-x4
855 dy1 = y(n1)-y4
856 dy2 = y(n2)-y4
857 dy3 = y(n3)-y4
858 dz1 = z(n1)-z4
859 dz2 = z(n2)-z4
860 dz3 = z(n3)-z4
861 
862 ! N4 lies above the plane of (N1,N2,N3) iff N3 lies above the plane of (N2,N1,N4) iff Det(N3-N4,N2-N4,N1-N4) = (N3-N4,N2-N4 X N1-N4) > 0.
863 call det(dx2,dy2,dz2,dx1,dy1,dz1,dx3,dy3,dz3,zz)
864 output = sup(zz,zero)
865 
866 ! Probe out
867 
868 
869 end subroutine stripack_swptst
870 
871 !----------------------------------------------------------------------
872 ! Subroutine: stripack_trfind
873 !> Locate a point relative to a triangulation
874 !----------------------------------------------------------------------
875 subroutine stripack_trfind(nst,p,n,x,y,z,list,lptr,lend,b1,b2,b3,i1,i2,i3)
876 
877 implicit none
878 
879 integer,intent(in) :: nst !< Index of a node at which TRFIND begins its search. Search time depends on the proximity of this node to P.
880 real(kind_real) :: p(3) !< X, Y, and Z coordinates (in that order) of the point P to be located.
881 integer,intent(in) :: n !< Number of nodes in the triangulation 3<=N.
882 real(kind_real),intent(in) :: x(n) !< X-coordinate of the triangulation nodes (unit vector).
883 real(kind_real),intent(in) :: y(n) !< Y-coordinate of the triangulation nodes (unit vector).
884 real(kind_real),intent(in) :: z(n) !< Z-coordinate of the triangulation nodes (unit vector).
885 integer,intent(in) :: list(*) !< Triangulation data structure created by TRMESH.
886 integer,intent(in) :: lptr(*) !< cf. list
887 integer,intent(in) :: lend(*) !< cf. list
888 real(kind_real),intent(out) :: b1 !< First unnormalized barycentric coordinate of the central projection of P onto the underlying planar triangle if P is in the convex hull of the nodes. These parameters are not altered if I1 = 0.
889 real(kind_real),intent(out) :: b2 !< Second unnormalized barycentric coordinate.
890 real(kind_real),intent(out) :: b3 !< Third unnormalized barycentric coordinate.
891 integer,intent(out) :: i1 !< First counterclockwise-ordered vertex index of a triangle containing P if P is contained in a triangle. If P is not in the convex hull of the nodes, I1 and I2 are the rightmost and leftmost (boundary) nodes that are visible from P, and I3 = 0. (If all boundary nodes are visible from P, then I1 and I2 coincide.) I1 = I2 = I3 = 0 if P and all of the nodes are coplanar (lie on a common great circle).
892 integer,intent(out) :: i2 !< Second counterclockwise-ordered vertex index.
893 integer,intent(out) :: i3 !< Third counterclockwise-ordered vertex index.
894 
895 ! Local variables
896 integer :: lp,n0,n1,n1s,n2,n2s,n3,n4,next,nf,nl
897 integer,save :: ix = 1,iy = 2,iz = 3
898 real(kind_real) :: output,eps,ptn1,ptn2,q(3),s12,tol,xp,yp,zp
899 
900 ! Set name
901 
902 
903 ! Probe in
904 
905 
906 ! Initialize variables.
907 xp = p(1)
908 yp = p(2)
909 zp = p(3)
910 n0 = nst
911 if ((n0<1).or.(n<n0)) call jrand(n,ix,iy,iz,n0)
912 
913 ! Compute the relative machine precision EPS and TOL.
914 eps = epsilon(eps)
915 tol = hundred*eps
916 
917 ! Set NF and NL to the first and last neighbors of N0, and initialize N1 = NF.
918 2 continue
919 lp = lend(n0)
920 nl = list(lp)
921 lp = lptr(lp)
922 nf = list(lp)
923 n1 = nf
924 
925 ! Find a pair of adjacent neighbors N1,N2 of N0 that define a wedge containing P: P LEFT N0->N1 and P RIGHT N0->N2.
926 if (0<nl) then
927  ! N0 is an interior node. Find N1.
928  3 continue
929  call det(x(n0),y(n0),z(n0),x(n1),y(n1),z(n1),xp,yp,zp,output)
930  if (inf(output,zero)) then
931  lp = lptr(lp)
932  n1 = list(lp)
933  if (n1 == nl) go to 6
934  go to 3
935  end if
936 else
937  ! N0 is a boundary node. Test for P exterior.
938  nl = -nl
939 
940  ! Is P to the right of the boundary edge N0->NF?
941  call det(x(n0),y(n0),z(n0),x(nf),y(nf),z(nf),xp,yp,zp,output)
942  if (inf(output,zero)) then
943  n1 = n0
944  n2 = nf
945  go to 9
946  end if
947 
948  ! Is P to the right of the boundary edge NL->N0?
949  call det(x(nl),y(nl),z(nl),x(n0),y(n0),z(n0),xp,yp,zp,output)
950  if (inf(output,zero)) then
951  n1 = nl
952  n2 = n0
953  go to 9
954  end if
955 end if
956 
957 ! P is to the left of arcs N0->N1 and NL->N0. Set N2 to the next neighbor of N0 (following N1).
958 4 continue
959 lp = lptr(lp)
960 n2 = abs(list(lp))
961 call det(x(n0),y(n0),z(n0),x(n2),y(n2),z(n2),xp,yp,zp,output)
962 if (inf(output,zero)) go to 7
963 n1 = n2
964 if (n1/=nl) go to 4
965 call det(x(n0),y(n0),z(n0),x(nf),y(nf),z(nf),xp,yp,zp,output)
966 if (inf(output,zero)) go to 6
967 
968 ! P is left of or on arcs N0->NB for all neighbors NB of N0. Test for P = +/-N0.
969 if (inf(abs(x(n0)*xp+y(n0)*yp+z(n0)*zp),one-four*eps)) then
970  ! All points are collinear iff P Left NB->N0 for all neighbors NB of N0. Search the neighbors of N0. Note: N1 = NL and LP points to NL.
971  do
972  call det(x(n1),y(n1),z(n1),x(n0),y(n0),z(n0),xp,yp,zp,output)
973  if (inf(output,zero)) exit
974  lp = lptr(lp)
975  n1 = abs(list(lp))
976  if (n1==nl) then
977  i1 = 0
978  i2 = 0
979  i3 = 0
980 
981  return
982  end if
983  end do
984 end if
985 
986 ! P is to the right of N1->N0, or P = +/-N0. Set N0 to N1 and start over.
987 n0 = n1
988 go to 2
989 
990 ! P is between arcs N0->N1 and N0->NF.
991 6 continue
992 n2 = nf
993 
994 ! P is contained in a wedge defined by geodesics N0-N1 and N0-N2, where N1 is adjacent to N2. Save N1 and N2 to test for cycling.
995 7 continue
996 n3 = n0
997 n1s = n1
998 n2s = n2
999 
1000 ! Top of edge-hopping loop:
1001 8 continue
1002 call det(x(n1),y(n1),z(n1),x(n2),y(n2),z(n2),xp,yp,zp,b3)
1003 if (inf(b3,zero)) then
1004  ! Set N4 to the first neighbor of N2 following N1 (the node opposite N2->N1) unless N1->N2 is a boundary arc.
1005  call lstptr(lend(n2),n1,list,lptr,lp)
1006  if (list(lp)<0) go to 9
1007  lp = lptr(lp)
1008  n4 = abs(list(lp))
1009 
1010  ! Define a new arc N1->N2 which intersects the geodesic N0-P.
1011  call det(x(n0),y(n0),z(n0),x(n4),y(n4),z(n4),xp,yp,zp,output)
1012  if (inf(output,zero)) then
1013  n3 = n2
1014  n2 = n4
1015  n1s = n1
1016  if ((n2/=n2s).and.(n2/=n0)) go to 8
1017  else
1018  n3 = n1
1019  n1 = n4
1020  n2s = n2
1021  if ((n1/=n1s).and.(n1/=n0)) go to 8
1022  end if
1023 
1024  ! The starting node N0 or edge N1-N2 was encountered again, implying a cycle (infinite loop). Restart with N0 randomly selected.
1025  call jrand(n,ix,iy,iz,n0)
1026  go to 2
1027 end if
1028 
1029 ! P is in (N1,N2,N3) unless N0, N1, N2, and P are collinear or P is close to -N0.
1030 if (supeq(b3,eps)) then
1031  ! B3 /= 0.
1032  call det(x(n2),y(n2),z(n2),x(n3),y(n3),z(n3),xp,yp,zp,b1)
1033  call det(x(n3),y(n3),z(n3),x(n1),y(n1),z(n1),xp,yp,zp,b2)
1034 
1035  ! Restart with N0 randomly selected.
1036  if (inf(b1,-tol).or.inf(b2,-tol)) then
1037  call jrand(n,ix,iy,iz,n0)
1038  go to 2
1039  end if
1040 else
1041  ! B3 = 0 and thus P lies on N1->N2. Compute B1 = Det(P,N2 X N1,N2) and B2 = Det(P,N1,N2 X N1).
1042  b3 = zero
1043  s12 = x(n1)*x(n2)+y(n1)*y(n2)+z(n1)*z(n2)
1044  ptn1 = xp*x(n1)+yp*y(n1)+zp*z(n1)
1045  ptn2 = xp*x(n2)+yp*y(n2)+zp*z(n2)
1046  b1 = ptn1-s12*ptn2
1047  b2 = ptn2-s12*ptn1
1048 
1049  ! Restart with N0 randomly selected.
1050  if (inf(b1,-tol).or.inf(b2,-tol)) then
1051  call jrand(n,ix,iy,iz,n0)
1052  go to 2
1053  end if
1054 end if
1055 
1056 ! P is in (N1,N2,N3).
1057 i1 = n1
1058 i2 = n2
1059 i3 = n3
1060 b1 = max(b1,zero)
1061 b2 = max(b2,zero)
1062 
1063 return
1064 
1065 ! P Right N1->N2, where N1->N2 is a boundary edge. Save N1 and N2, and set NL = 0 to indicate that NL has not yet been found.
1066 9 continue
1067 n1s = n1
1068 n2s = n2
1069 nl = 0
1070 
1071 ! Counterclockwise Boundary Traversal:
1072 10 continue
1073 lp = lend(n2)
1074 lp = lptr(lp)
1075 next = list(lp)
1076 call det(x(n2),y(n2),z(n2),x(next),y(next),z(next),xp,yp,zp,output)
1077 if (supeq(output,zero)) then
1078  ! N2 is the rightmost visible node if P Forward N2->N1 or NEXT Forward N2->N1. Set Q to (N2 X N1) X N2.
1079  s12 = x(n1)*x(n2)+y(n1)*y(n2)+z(n1)*z(n2)
1080  q(1) = x(n1)-s12*x(n2)
1081  q(2) = y(n1)-s12*y(n2)
1082  q(3) = z(n1)-s12*z(n2)
1083  if (sup(xp*q(1)+yp*q(2)+zp*q(3),zero).or.sup(x(next)*q(1)+y(next)*q(2)+z(next)*q(3),zero)) go to 11
1084 
1085  ! N1, N2, NEXT, and P are nearly collinear, and N2 is the leftmost visible node.
1086  nl = n2
1087 end if
1088 
1089 ! Bottom of counterclockwise loop:
1090 n1 = n2
1091 n2 = next
1092 if (n2/=n1s) go to 10
1093 
1094 ! All boundary nodes are visible from P.
1095 i1 = n1s
1096 i2 = n1s
1097 i3 = 0
1098 
1099 return
1100 
1101 ! N2 is the rightmost visible node.
1102 11 continue
1103 nf = n2
1104 if (nl==0) then
1105  ! Restore initial values of N1 and N2, and begin the search for the leftmost visible node.
1106  n2 = n2s
1107  n1 = n1s
1108 
1109  ! Clockwise Boundary Traversal:
1110  12 continue
1111  lp = lend(n1)
1112  next = -list(lp)
1113  call det(x(next),y(next),z(next),x(n1),y(n1),z(n1),xp,yp,zp,output)
1114  if (supeq(output,zero)) then
1115  ! N1 is the leftmost visible node if P or NEXT is forward of N1->N2. Compute Q = N1 X (N2 X N1).
1116  s12 = x(n1)*x(n2)+y(n1)*y(n2)+z(n1)*z(n2)
1117  q(1) = x(n2)-s12*x(n1)
1118  q(2) = y(n2)-s12*y(n1)
1119  q(3) = z(n2)-s12*z(n1)
1120  if (sup(xp*q(1)+yp*q(2)+zp*q(3),zero).or.sup(x(next)*q(1)+y(next)*q(2)+z(next)*q(3),zero)) go to 13
1121 
1122  ! P, NEXT, N1, and N2 are nearly collinear and N1 is the rightmost visible node.
1123  nf = n1
1124  end if
1125 
1126  ! Bottom of clockwise loop:
1127  n2 = n1
1128  n1 = next
1129  if (n1/=n1s) go to 12
1130 
1131  ! All boundary nodes are visible from P.
1132  i1 = n1
1133  i2 = n1
1134  i3 = 0
1135 
1136  return
1137 
1138  ! N1 is the leftmost visible node.
1139  13 continue
1140  nl = n1
1141 end if
1142 
1143 ! NF and NL have been found.
1144 i1 = nf
1145 i2 = nl
1146 i3 = 0
1147 
1148 ! Probe out
1149 
1150 
1151 end subroutine stripack_trfind
1152 
1153 !----------------------------------------------------------------------
1154 ! Subroutine: stripack_trmesh
1155 !> Create a Delaunay triangulation on the unit sphere
1156 !----------------------------------------------------------------------
1157 subroutine stripack_trmesh(mpl,n,x,y,z,list,lptr,lend,lnew,near,next,dist)
1158 
1159 implicit none
1160 
1161 type(mpl_type),intent(inout) :: mpl !< MPI data
1162 integer,intent(in) :: n !< Number of nodes in the triangulation. 3<=N.
1163 real(kind_real),intent(in) :: x(n) !< X-coordinate of distinct nodes. (X(K),Y(K), Z(K)) is referred to as node K, and K is referred to as a nodal index. It is required that X(K)**2+Y(K)**2+Z(K)**2 = 1 for all K. The first three nodes must not be colinear (lie on a common great circle).
1164 real(kind_real),intent(in) :: y(n) !< Y-coordinate of distinct nodes.
1165 real(kind_real),intent(in) :: z(n) !< Z-coordinate of distinct nodes.
1166 integer,intent(out) :: list(6*(n-2)) !< Nodal indexes which, along with LPTR, LEND, and LNEW, define the triangulation as a set of N adjacency lists; counterclockwise-ordered sequences of neighboring nodes such that the first and last neighbors of a boundary node are boundary nodes (the first neighbor of an interior node is arbitrary). In order to distinguish between interior and boundary nodes, the last neighbor of each boundary node is represented by the negative of its index.
1167 integer,intent(out) :: lptr(6*(n-2)) !< Set of pointers (LIST indexes) in one-to-one correspondence with the elements of LIST. LIST(LPTR(I)) indexes the node which follows LIST(I) in cyclical counterclockwise order (the first neighbor follows the last neighbor).
1168 integer,intent(out) :: lend(n) !< Pointers to adjacency lists. LEND(K) points to the last neighbor of node K. LIST(LEND(K))<0 if and only if K is a boundary node.
1169 integer,intent(out) :: lnew !< Pointer to the first empty location in LIST and LPTR (list length plus one).
1170 integer,intent(inout) :: near(n) !< Workspace used to efficiently determine the nearest triangulation node to each unprocessed node for use by ADDNOD.
1171 integer,intent(inout) :: next(n) !< Workspace used to efficiently determine the nearest triangulation node to each unprocessed node for use by ADDNOD.
1172 real(kind_real),intent(inout) :: dist(n) !< Workspace used to efficiently determine the nearest triangulation node to each unprocessed node for use by ADDNOD.
1173 
1174 ! Local variables
1175 integer :: i,i0,j,k,lp,lpl,nexti,nn
1176 real(kind_real) :: d,d1,d2,d3
1177 logical :: output_1,output_2
1178 
1179 ! Set name
1180 
1181 
1182 ! Probe in
1183 
1184 
1185 ! Initialization
1186 nn = n
1187 if (nn<3) call mpl%abort('stripack_trmesh','N < 3')
1188 list = 0
1189 lptr = 0
1190 lend = 0
1191 
1192 ! Store the first triangle in the linked list.
1193 call left(x(1),y(1),z(1),x(2),y(2),z(2),x(3),y(3),z(3),output_1)
1194 call left(x(2),y(2),z(2),x(1),y(1),z(1),x(3),y(3),z(3),output_2)
1195 
1196 if (.not.output_1) then
1197  ! The first triangle is (3,2,1) = (2,1,3) = (1,3,2).
1198  list(1) = 3
1199  lptr(1) = 2
1200  list(2) = -2
1201  lptr(2) = 1
1202  lend(1) = 2
1203 
1204  list(3) = 1
1205  lptr(3) = 4
1206  list(4) = -3
1207  lptr(4) = 3
1208  lend(2) = 4
1209 
1210  list(5) = 2
1211  lptr(5) = 6
1212  list(6) = -1
1213  lptr(6) = 5
1214  lend(3) = 6
1215 else if (.not.output_2) then
1216  ! The first triangle is (1,2,3): 3 Strictly Left 1->2, i.e., node 3 lies in the left hemisphere defined by arc 1->2.
1217  list(1) = 2
1218  lptr(1) = 2
1219  list(2) = -3
1220  lptr(2) = 1
1221  lend(1) = 2
1222 
1223  list(3) = 3
1224  lptr(3) = 4
1225  list(4) = -1
1226  lptr(4) = 3
1227  lend(2) = 4
1228 
1229  list(5) = 1
1230  lptr(5) = 6
1231  list(6) = -2
1232  lptr(6) = 5
1233  lend(3) = 6
1234 else
1235  ! The first three nodes are colinear.
1236  call mpl%abort('stripack_trmesh','The first three nodes are colinear')
1237 end if
1238 
1239 ! Initialize LNEW and test for N = 3.
1240 lnew = 7
1241 if (nn==3) then
1242 
1243  return
1244 end if
1245 
1246 ! A nearest-node data structure (NEAR, NEXT, and DIST) is used to obtain an expected-time (N*log(N)) incremental algorithm by enabling constant search time for locating each new node in the triangulation. For each unprocessed node K, NEAR(K) is the index of the triangulation node closest to K (used as the starting point for the search in Subroutine TRFIND) and DIST(K) is an increasing function of the arc length (angular distance) between nodes K and NEAR(K): -cos(a) for arc length a. Since it is necessary to efficiently find the subset of unprocessed nodes associated with each triangulation node J (those that have J as their NEAR entries), the subsets are stored in NEAR and NEXT as follows: for each node J in the triangulation, I = NEAR(J) is the first unprocessed node in J's set (with I = 0 if the set is empty), L = NEXT(I) (if 0<I) is the second, NEXT(L) (if 0<L) is the third, etc. The nodes in each set are initially ordered by increasing indexes (which maximizes efficiency) but that ordering is not maintained as the data structure is updated.
1247 
1248 ! Initialize the data structure for the single triangle.
1249 near(1) = 0
1250 near(2) = 0
1251 near(3) = 0
1252 do k=nn,4,-1
1253  d1 = -(x(k)*x(1)+y(k)*y(1)+z(k)*z(1))
1254  d2 = -(x(k)*x(2)+y(k)*y(2)+z(k)*z(2))
1255  d3 = -(x(k)*x(3)+y(k)*y(3)+z(k)*z(3))
1256  if (inf(d1,d2).and.inf(d1,d3)) then
1257  near(k) = 1
1258  dist(k) = d1
1259  next(k) = near(1)
1260  near(1) = k
1261  else if (inf(d2,d1).and.inf(d2,d3)) then
1262  near(k) = 2
1263  dist(k) = d2
1264  next(k) = near(2)
1265  near(2) = k
1266  else
1267  near(k) = 3
1268  dist(k) = d3
1269  next(k) = near(3)
1270  near(3) = k
1271  end if
1272 end do
1273 
1274 ! Add the remaining nodes.
1275 do k=4,nn
1276  call addnod(mpl,near(k),k,x,y,z,list,lptr,lend,lnew)
1277 
1278  ! Remove K from the set of unprocessed nodes associated with NEAR(K).
1279  i = near(k)
1280  i0 = i
1281  if (near(i)==k) then
1282  near(i) = next(k)
1283  else
1284  i = near(i)
1285  do
1286  i0 = i
1287  i = next(i0)
1288  if (i==k) exit
1289  end do
1290  next(i0) = next(k)
1291  end if
1292  near(k) = 0
1293 
1294  ! Loop on neighbors J of node K.
1295  lpl = lend(k)
1296  lp = lpl
1297  3 continue
1298  lp = lptr(lp)
1299  j = abs(list(lp))
1300 
1301  ! Loop on elements I in the sequence of unprocessed nodes associated with J: K is a candidate for replacing J as the nearest triangulation node to I. The next value of I in the sequence, NEXT(I), must be saved before I is moved because it is altered by adding I to K's set.
1302  i = near(j)
1303  do
1304  if (i==0) exit
1305  nexti = next(i)
1306 
1307  ! Test for the distance from I to K less than the distance from I to J. Indistinguishability threshold for cross-platform reproducibility
1308  d =-(x(i)*x(k)+y(i)*y(k)+z(i)*z(k))
1309  if (inf(d,dist(i))) then
1310  ! Replace J by K as the nearest triangulation node to I: update NEAR(I) and DIST(I), and remove I from J's set of unprocessed nodes and add it to K's set.
1311  near(i) = k
1312  dist(i) = d
1313  if (i==near(j)) then
1314  near(j) = nexti
1315  else
1316  next(i0) = nexti
1317  end if
1318  next(i) = near(k)
1319  near(k) = i
1320  else
1321  i0 = i
1322  end if
1323  i = nexti
1324  end do
1325 
1326  ! Bottom of loop on neighbors J.
1327  if (lp/=lpl) go to 3
1328 end do
1329 
1330 ! Probe out
1331 
1332 
1333 end subroutine stripack_trmesh
1334 
1335 end module tools_stripack
Subroutines/functions list.
Definition: tools_const.F90:31
real(kind_real), parameter, public one
One.
Definition: tools_const.F90:42
real(kind_real), parameter, public pi
Pi.
Definition: tools_const.F90:53
real(kind_real), parameter, public zero
Zero.
Definition: tools_const.F90:37
real(kind_real), parameter, public hundred
Hundred.
Definition: tools_const.F90:51
real(kind_real), parameter, public four
Four.
Definition: tools_const.F90:45
Kinds definition.
Definition: tools_kinds.F90:9
integer, parameter, public kind_real
Real kind alias for the whole code.
Definition: tools_kinds.F90:31
Generic ranks, dimensions and types.
Definition: tools_repro.F90:42
Subroutines/functions list.
subroutine stripack_trmesh(mpl, n, x, y, z, list, lptr, lend, lnew, near, next, dist)
Create a Delaunay triangulation on the unit sphere.
subroutine stripack_bdyadd(kk, i1, i2, list, lptr, lend, lnew)
Add a boundary node to a triangulation.
subroutine stripack_swap(in1, in2, io1, io2, list, lptr, lend, lp21)
Replace the diagonal arc of a quadrilateral with the other diagonal.
subroutine stripack_insert(k, lp, list, lptr, lnew)
Insert K as a neighbor of N1.
subroutine stripack_jrand(n, ix, iy, iz, output)
Return a random integer between 1 and N.
subroutine stripack_covsph(kk, n0, list, lptr, lend, lnew)
Connect an exterior node to boundary nodes, covering the sphere.
subroutine stripack_area(v1, v2, v3, area)
Compute the area of a spherical triangle.
subroutine stripack_bnodes(n, list, lptr, lend, nodes, nb)
Return the boundary nodes of a triangulation.
subroutine stripack_intadd(kk, i1, i2, i3, list, lptr, lend, lnew)
Add an interior node to a triangulation.
subroutine stripack_trfind(nst, p, n, x, y, z, list, lptr, lend, b1, b2, b3, i1, i2, i3)
Locate a point relative to a triangulation.
subroutine stripack_addnod(mpl, nst, k, x, y, z, list, lptr, lend, lnew)
Add a node to a triangulation.
subroutine stripack_swptst(n1, n2, n3, n4, x, y, z, output)
Decide whether to replace a diagonal arc by the other.
subroutine stripack_det(x1, y1, z1, x2, y2, z2, x0, y0, z0, output)
Compute 3D determinant.
subroutine stripack_left(x1, y1, z1, x2, y2, z2, x0, y0, z0, output)
Determine whether a node is to the left of a plane through the origin.
subroutine stripack_lstptr(lpl, nb, list, lptr, output)
Return the index of NB in the adjacency list.
Generic ranks, dimensions and types.
Definition: type_mpl.F90:42