@@ -395,6 +395,87 @@ def from_file( cls, filename ):
395
395
396
396
@classmethod
397
397
def join (cls , patches , connectivity , name ):
398
+ """
399
+ creates a multipatch domain by joining two or more patches in 2D or 3D
400
+
401
+ Parameters
402
+ ----------
403
+ patches : list
404
+ list of patches
405
+
406
+ connectivity : list
407
+ list of interfaces, identified by a tuple of 2 boundaries and an orientation: (bound_minus, bound_plus, ornt)
408
+ where
409
+ - each boundary is identified by a tuple of 3 integers: (patch, axis, ext)
410
+ with patches given as objects (or by their indices in the patches list)
411
+ and
412
+ - In 2D, ornt is an integer that can take the value of 1 or -1
413
+ - In 3D, ornt is a tuple of 3 integers that can take the value of 1 or -1
414
+ (see below for more details)
415
+
416
+
417
+ name : str
418
+ name of the domain
419
+
420
+ Returns
421
+ -------
422
+ domain : Domain
423
+ multipatch domain
424
+
425
+ Notes
426
+ -----
427
+ The orientations are specified in the same manner as in GeoPDES, see e.g.
428
+ <https://github.com/rafavzqz/geopdes/blob/master/geopdes/doc/geo_specs_mp_v21.txt#L193-L237>
429
+ and
430
+ T. Dokken, E. Quak, V. Skytt. Requirements from Isogeometric Analysis for changes in product design ontologies, 2010.
431
+
432
+ Example
433
+ -------
434
+ # list of patches (mapped domains)
435
+ Omega_0 = F0(A)
436
+ Omega_1 = F1(A)
437
+ Omega_2 = F2(A)
438
+ Omega_3 = F3(A)
439
+
440
+ patches = [Omega_0, Omega_1, Omega_2, Omega_3]
441
+
442
+ # integers representing the axes
443
+ axis_0 = 0
444
+ axis_1 = 1
445
+ axis_2 = 2
446
+
447
+ # integers representing the extremities: left (-1) or right (+1)
448
+ ext_0 = -1
449
+ ext_1 = +1
450
+
451
+ # A connectivity list in 2D
452
+ connectivity = [((Omega_0, axis_0, ext_0), (Omega_1, axis_0, ext_1), 1),
453
+ ((Omega_1, axis_1, ext_0), (Omega_3, axis_1, ext_1), -1),
454
+ ((Omega_0, axis_1, ext_0), (Omega_2, axis_1, ext_1), 1),
455
+ ((Omega_2, axis_0, ext_0), (Omega_3, axis_0, ext_1), -1)]
456
+
457
+ # alternative option (passing interface patches by their indices in the patches list):
458
+ connectivity = [((0, axis_0, ext_0), (1, axis_0, ext_1), 1),
459
+ ((1, axis_1, ext_0), (3, axis_1, ext_1), -1),
460
+ ((0, axis_1, ext_0), (2, axis_1, ext_1), 1),
461
+ ((2, axis_0, ext_0), (3, axis_0, ext_1), -1)]
462
+
463
+ # A connectivity list in 3D
464
+ connectivity = [((Omega_0, axis_0, ext_1), (Omega_1, axis_0, ext_0), ( 1, 1, 1)),
465
+ ((Omega_0, axis_1, ext_1), (Omega_2, axis_1, ext_0), ( 1, -1, 1)),
466
+ ((Omega_1, axis_1, ext_1), (Omega_3, axis_1, ext_0), (-1, 1, -1)),
467
+ ((Omega_2, axis_0, ext_1), (Omega_3, axis_0, ext_0), (-1, 1, 1))]
468
+
469
+ # alternative option (passing interface patches by their indices in the patches list):
470
+ connectivity = [((0, axis_0, ext_1), (1, axis_0, ext_0), ( 1, 1, 1)),
471
+ ((0, axis_1, ext_1), (2, axis_1, ext_0), ( 1, -1, 1)),
472
+ ((1, axis_1, ext_1), (3, axis_1, ext_0), (-1, 1, -1)),
473
+ ((2, axis_0, ext_1), (3, axis_0, ext_0), (-1, 1, 1))]
474
+
475
+ # the multi-patch domain
476
+ Omega = Domain.join(patches=patches, connectivity=connectivity, name='Omega')
477
+
478
+ """
398
479
399
480
assert isinstance (patches , (tuple , list ))
400
481
assert isinstance (connectivity , (tuple , list ))
@@ -403,13 +484,21 @@ def join(cls, patches, connectivity, name):
403
484
assert all (p .dim == patches [0 ].dim for p in patches )
404
485
dim = int (patches [0 ].dim )
405
486
487
+ patch_given_by_indices = (len (connectivity ) > 0 and isinstance (connectivity [0 ][0 ][0 ], int ))
488
+
406
489
from sympde .topology .mapping import MultiPatchMapping
407
490
# ... connectivity
408
491
interfaces = {}
409
492
boundaries = []
410
493
for cn in connectivity :
411
- bnd_minus = patches [cn [0 ][0 ]].get_boundary (axis = cn [0 ][1 ], ext = cn [0 ][2 ])
412
- bnd_plus = patches [cn [1 ][0 ]].get_boundary (axis = cn [1 ][1 ], ext = cn [1 ][2 ])
494
+ if patch_given_by_indices :
495
+ patch_minus = patches [cn [0 ][0 ]]
496
+ patch_plus = patches [cn [1 ][0 ]]
497
+ else :
498
+ patch_minus = cn [0 ][0 ]
499
+ patch_plus = cn [1 ][0 ]
500
+ bnd_minus = patch_minus .get_boundary (axis = cn [0 ][1 ], ext = cn [0 ][2 ])
501
+ bnd_plus = patch_plus .get_boundary (axis = cn [1 ][1 ], ext = cn [1 ][2 ])
413
502
if dim == 1 : ornt = None
414
503
elif dim == 2 : ornt = 1 if len (cn ) == 2 else cn [2 ]
415
504
elif dim == 3 :
0 commit comments