Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/fem/src/LumpingUtils.F90
5218 views
1
!/*****************************************************************************/
2
! *
3
! * Elmer, A Finite Element Software for Multiphysical Problems
4
! *
5
! * Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland
6
! *
7
! * This library is free software; you can redistribute it and/or
8
! * modify it under the terms of the GNU Lesser General Public
9
! * License as published by the Free Software Foundation; either
10
! * version 2.1 of the License, or (at your option) any later version.
11
! *
12
! * This library is distributed in the hope that it will be useful,
13
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
14
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
15
! * Lesser General Public License for more details.
16
! *
17
! * You should have received a copy of the GNU Lesser General Public
18
! * License along with this library (in file ../LGPL-2.1); if not, write
19
! * to the Free Software Foundation, Inc., 51 Franklin Street,
20
! * Fifth Floor, Boston, MA 02110-1301 USA
21
! *
22
! *****************************************************************************/
23
!
24
!/******************************************************************************
25
! *
26
! * Generic utilities related to components
27
! *
28
! ******************************************************************************
29
! *
30
! * Authors: Juha Ruokolainen
31
! * Email: [email protected]
32
! * Web: http://www.csc.fi/elmer
33
! * Address: CSC - IT Center for Science Ltd.
34
! * Keilaranta 14
35
! * 02101 Espoo, Finland
36
! *
37
! * Original Date: 28 Sep 1998
38
! *
39
! *****************************************************************************/
40
41
!> Generic utilities related to components
42
!------------------------------------------------------------------------------
43
44
!> \ingroup ElmerLib
45
!> \{
46
47
48
MODULE LumpingUtils
49
50
USE Lists
51
USE ElementUtils
52
USE ElementDescription
53
USE ParallelUtils
54
55
IMPLICIT NONE
56
57
58
CONTAINS
59
60
61
!------------------------------------------------------------------------------
62
!> Compute reduction operators for a given component with given nodal vector field.
63
!> Force is simple sum of nodal forces
64
!> Moment is moment of nodal forces about a given center point
65
!> Torque is moment of nodal forces about a given rotational axis
66
!> If the given field is a elemental (DG) field it may be reduced by
67
!> optional SetPerm reordering for minimal discontinuous set.
68
!------------------------------------------------------------------------------
69
SUBROUTINE ComponentNodalForceReduction(Model, Mesh, CompParams, NF, &
70
Force, Moment, Torque, SetPerm )
71
!------------------------------------------------------------------------------
72
TYPE(Model_t) :: Model
73
TYPE(Mesh_t), POINTER :: Mesh
74
TYPE(ValueList_t), POINTER :: CompParams
75
TYPE(Variable_t), POINTER :: NF
76
REAL(KIND=dp), OPTIONAL :: Moment(3), Force(3), Torque
77
INTEGER, POINTER, OPTIONAL :: SetPerm(:)
78
!------------------------------------------------------------------------------
79
! Local variables
80
!------------------------------------------------------------------------------
81
TYPE(Element_t), POINTER :: Element
82
TYPE(Nodes_t) :: Nodes
83
LOGICAL, ALLOCATABLE :: VisitedNode(:)
84
REAL(KIND=dp) :: Origin(3), Axis(3), P(3), F(3), v1(3), v2(3), &
85
RotorRadius, rad, minrad, maxrad, eps
86
REAL(KIND=dp), POINTER :: Pwrk(:,:)
87
INTEGER :: t, i, j, k, n, dofs, globalnode, AirBody
88
LOGICAL :: ElementalVar, Found, NeedLocation
89
INTEGER, POINTER :: MasterEntities(:),NodeIndexes(:),DofIndexes(:)
90
LOGICAL :: VisitNodeOnlyOnce
91
INTEGER :: FirstElem, LastElem
92
LOGICAL :: BcMode, BulkMode, RotorMode, isParallel
93
CHARACTER(*), PARAMETER :: Caller = 'ComponentNodalForceReduction'
94
95
CALL Info(Caller,'Performing reduction for component: '&
96
//TRIM(ListGetString(CompParams,'Name')),Level=10)
97
98
IF(.NOT. (PRESENT(Torque) .OR. PRESENT(Moment) .OR. PRESENT(Force) ) ) THEN
99
CALL Warn(Caller,'Nothing to compute!')
100
RETURN
101
END IF
102
103
IF( PRESENT(Torque)) Torque = 0.0_dp
104
IF( PRESENT(Moment)) Moment = 0.0_dp
105
IF( PRESENT(Force)) Force = 0.0_dp
106
107
isParallel = CurrentModel % Solver % Parallel
108
109
110
eps = 1.0e-6
111
BcMode = .FALSE.
112
BulkMode = .FALSE.
113
RotorMode = ListGetLogical( CompParams,'Rotor Mode',Found )
114
IF( RotorMode ) THEN
115
RotorRadius = ListGetConstReal( CurrentModel % Simulation,'Rotor Radius',Found )
116
IF(.NOT. Found ) THEN
117
CALL Fatal(Caller,'"Rotor Mode" requires "Rotor Radius"')
118
END IF
119
ELSE
120
MasterEntities => ListGetIntegerArray( CompParams,'Master Bodies',BulkMode )
121
IF( .NOT. BulkMode ) THEN
122
MasterEntities => ListGetIntegerArray( CompParams,'Master Boundaries', BCMode)
123
END IF
124
IF(.NOT. (BulkMode .OR. BCMode ) ) THEN
125
CALL Warn(Caller,'> Master Bodies < or > Master Boundaries < not given')
126
RETURN
127
END IF
128
END IF
129
130
NeedLocation = PRESENT( Moment ) .OR. PRESENT( Torque )
131
132
! User may specific origin and axis for torque computation
133
! By default (0,0,0) is the origin, and (0,0,1) the axis.
134
Pwrk => ListGetConstRealArray( CompParams,'Torque Origin',Found )
135
IF( Found ) THEN
136
IF( SIZE(Pwrk,1) /= 3 .OR. SIZE(Pwrk,2) /= 1 ) THEN
137
CALL Fatal(Caller,'Size of > Torque Origin < should be 3!')
138
END IF
139
Origin = Pwrk(1:3,1)
140
ELSE
141
Origin = 0.0_dp
142
END IF
143
Pwrk => ListGetConstRealArray( CompParams,'Torque Axis',Found )
144
IF( Found ) THEN
145
IF( SIZE(Pwrk,1) /= 3 .OR. SIZE(Pwrk,2) /= 1 ) THEN
146
CALL Fatal(Caller,'Size of > Torque Axis < should be 3!')
147
END IF
148
Axis = Pwrk(1:3,1)
149
! Normalize axis is it should just be used for the direction
150
Axis = Axis / SQRT( SUM( Axis*Axis ) )
151
ELSE
152
Axis = 0.0_dp
153
Axis(3) = 1.0_dp
154
END IF
155
156
ElementalVar = ( NF % TYPE == Variable_on_nodes_on_elements )
157
IF( PRESENT( SetPerm ) .AND. .NOT. ElementalVar ) THEN
158
CALL Fatal(Caller,'SetPerm is usable only for elemental fields')
159
END IF
160
161
dofs = NF % Dofs
162
IF( dofs == 2 ) F(3) = 0.0_dp
163
164
! For nodal field compute only once each node
165
! For DG field each node is visited only once by construction
166
VisitNodeOnlyOnce = .NOT. ElementalVar .OR. PRESENT(SetPerm)
167
IF( VisitNodeOnlyOnce ) THEN
168
IF( PRESENT( SetPerm ) ) THEN
169
n = MAXVAL( SetPerm )
170
ELSE
171
n = Mesh % NumberOfNodes
172
END IF
173
ALLOCATE(VisitedNode( n ) )
174
VisitedNode = .FALSE.
175
END IF
176
177
IF( BcMode ) THEN
178
FirstElem = Mesh % NumberOfBulkElements + 1
179
LastElem = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
180
ELSE
181
FirstElem = 1
182
LastElem = Mesh % NumberOfBulkElements
183
END IF
184
185
! This is a special reduction that only applies to rotors that are surrounded by airgap.
186
! Very special operator for electrical machines that removes the bookkeeping of the bodies
187
! that constitute the rotor.
188
AirBody = 0
189
IF(RotorMode ) THEN
190
AirBody = ListGetInteger( CompParams,'Air Body',Found )
191
IF(AirBody == 0) THEN
192
DO t=FirstElem,LastElem
193
Element => Mesh % Elements(t)
194
n = Element % TYPE % NumberOfNodes
195
CALL CopyElementNodesFromMesh( Nodes, Mesh, n, Element % NodeIndexes)
196
DO i=1,n
197
rad = SQRT(Nodes % x(i)**2 + Nodes % y(i)**2)
198
IF(i==1) THEN
199
minrad = rad
200
maxrad = rad
201
ELSE
202
minrad = MIN(minrad,rad)
203
maxrad = MAX(maxrad,rad)
204
END IF
205
END DO
206
207
! The body is defined by an element that is at and inside the rotor radius.
208
IF(ABS(maxrad-RotorRadius) < eps .AND. minrad < RotorRadius*(1-eps) ) THEN
209
AirBody = Element % BodyId
210
EXIT
211
END IF
212
END DO
213
AirBody = ParallelReduction(AirBody,2)
214
CALL Info(Caller,'Airgap inner body determined to be: '//I2S(AirBody),Level=12)
215
IF(AirBody==0) THEN
216
CALL Fatal(Caller,'Could not define airgap inner body!')
217
ELSE
218
CALL ListAddInteger(CompParams,'Air Body',AirBody)
219
END IF
220
END IF
221
END IF
222
223
224
DO t=FirstElem,LastElem
225
Element => Mesh % Elements(t)
226
227
IF( BcMode ) THEN
228
IF( ALL( MasterEntities /= Element % BoundaryInfo % Constraint ) ) CYCLE
229
ELSE IF( BulkMode ) THEN
230
IF( ALL( MasterEntities /= Element % BodyId ) ) CYCLE
231
ELSE IF( RotorMode ) THEN
232
IF( Element % BodyId == AirBody ) CYCLE
233
CALL CopyElementNodesFromMesh( Nodes, Mesh, n, Element % NodeIndexes)
234
rad = SQRT((SUM(Nodes % x(1:n))/n)**2 + (SUM(Nodes % y(1:n))/n)**2)
235
IF(rad > RotorRadius ) CYCLE
236
END IF
237
238
n = Element % TYPE % NumberOfNodes
239
NodeIndexes => Element % NodeIndexes
240
IF( ElementalVar ) THEN
241
DofIndexes => Element % DGIndexes
242
ELSE
243
DofIndexes => NodeIndexes
244
END IF
245
246
DO i=1,n
247
j = DofIndexes(i)
248
k = NF % Perm(j)
249
IF( k == 0 ) CYCLE
250
251
IF( VisitNodeOnlyOnce ) THEN
252
IF( PRESENT( SetPerm ) ) j = SetPerm(j)
253
IF( VisitedNode(j) ) CYCLE
254
VisitedNode(j) = .TRUE.
255
END IF
256
257
globalnode = NodeIndexes(i)
258
259
! Only compute the parallel reduction once
260
IF( isParallel ) THEN
261
IF( Element % PartIndex /= ParEnv % MyPe ) CYCLE
262
263
! This is (probably) not correct, the "nodal forces"-array is partial and should be summed --> comment out.
264
! IF( VisitNodeOnlyOnce ) THEN
265
! IF( Mesh % ParallelInfo % NeighbourList(globalnode) % Neighbours(1) /= ParEnv % MyPE ) CYCLE
266
! END IF
267
END IF
268
269
F(1) = NF % Values( dofs*(k-1) + 1)
270
F(2) = NF % Values( dofs*(k-1) + 2)
271
IF( dofs == 3 ) THEN
272
F(3) = NF % Values( dofs*(k-1) + 3)
273
END IF
274
275
IF( PRESENT( Force ) ) THEN
276
! calculate simple sum
277
Force = Force + F
278
END IF
279
280
IF( NeedLocation ) THEN
281
P(1) = Mesh % Nodes % x(globalnode)
282
P(2) = Mesh % Nodes % y(globalnode)
283
P(3) = Mesh % Nodes % z(globalnode)
284
285
v1 = P - Origin
286
287
! Calculate moment
288
IF( PRESENT( Moment ) ) THEN
289
Moment = Moment + CrossProduct(v1,F)
290
END IF
291
292
! Calculate torque around an axis
293
IF( PRESENT( Torque ) ) THEN
294
v1 = v1 - SUM(Axis*v1)*Axis
295
v2 = CrossProduct(v1,F)
296
Torque = Torque + SUM(Axis*v2)
297
END IF
298
END IF
299
300
END DO
301
END DO
302
303
IF( isParallel ) THEN
304
IF( PRESENT( Force ) ) THEN
305
DO i=1,3
306
Force(i) = ParallelReduction(Force(i))
307
END DO
308
END IF
309
310
IF( PRESENT( Moment ) ) THEN
311
DO i=1,3
312
Moment(i) = ParallelReduction(Moment(i))
313
END DO
314
END IF
315
316
IF( PRESENT( Torque ) ) THEN
317
Torque = ParallelReduction(Torque)
318
END IF
319
END IF
320
321
!------------------------------------------------------------------------------
322
END SUBROUTINE ComponentNodalForceReduction
323
!------------------------------------------------------------------------------
324
325
326
!------------------------------------------------------------------------------
327
!> Perform reduction from distributed data to components.
328
!> These are similar operations as the stastistical operations in SaveScalars.
329
!------------------------------------------------------------------------------
330
FUNCTION ComponentNodalReduction(Model, Mesh, CompParams, Var, OperName ) RESULT ( OperX )
331
!------------------------------------------------------------------------------
332
TYPE(Model_t) :: Model
333
TYPE(Mesh_t), POINTER :: Mesh
334
TYPE(ValueList_t), POINTER :: CompParams
335
TYPE(Variable_t), POINTER :: Var
336
REAL(KIND=dp) :: OperX
337
CHARACTER(LEN=*) :: OperName
338
!------------------------------------------------------------------------------
339
! Local variables
340
!------------------------------------------------------------------------------
341
TYPE(Element_t), POINTER :: Element
342
LOGICAL, ALLOCATABLE :: VisitedNode(:)
343
INTEGER :: t, i, j, k, l, n, NoDofs, globalnode, sumi
344
REAL(KIND=dp) :: X, Minimum, Maximum, AbsMinimum, AbsMaximum, SumX, SumXX, SumAbsX
345
LOGICAL :: ElementalVar, Found
346
INTEGER, POINTER :: MasterEntities(:),NodeIndexes(:),DofIndexes(:)
347
LOGICAL :: VisitNodeOnlyOnce, Initialized
348
INTEGER :: FirstElem, LastElem
349
LOGICAL :: BcMode
350
351
352
CALL Info('ComponentNodalReduction','Performing reduction for component: '&
353
//TRIM(ListGetString(CompParams,'Name')),Level=10)
354
355
OperX = 0.0_dp
356
357
BcMode = .FALSE.
358
MasterEntities => ListGetIntegerArray( CompParams,'Master Bodies',Found )
359
IF( .NOT. Found ) THEN
360
MasterEntities => ListGetIntegerArray( CompParams,'Master Boundaries',Found )
361
BcMode = .TRUE.
362
END IF
363
364
IF(.NOT. Found ) THEN
365
CALL Warn('ComponentNodalReduction',&
366
'> Master Bodies < or > Master Boundaries < not given')
367
RETURN
368
END IF
369
370
IF( BcMode ) THEN
371
FirstElem = Mesh % NumberOfBulkElements + 1
372
LastElem = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
373
ELSE
374
FirstElem = 1
375
LastElem = Mesh % NumberOfBulkElements
376
END IF
377
378
ElementalVar = ( Var % TYPE == Variable_on_nodes_on_elements )
379
NoDofs = Var % Dofs
380
Initialized = .FALSE.
381
382
! For nodal field compute only once each node
383
! For DG field each node is visited only once by construction
384
VisitNodeOnlyOnce = .NOT. ElementalVar
385
IF( VisitNodeOnlyOnce ) THEN
386
n = Mesh % NumberOfNodes
387
ALLOCATE(VisitedNode( n ) )
388
VisitedNode = .FALSE.
389
END IF
390
391
sumi = 0
392
sumx = 0.0_dp
393
sumxx = 0.0_dp
394
sumabsx = 0.0_dp
395
Maximum = 0.0_dp
396
Minimum = 0.0_dp
397
AbsMaximum = 0.0_dp
398
AbsMinimum = 0.0_dp
399
400
DO t=FirstElem,LastElem
401
Element => Mesh % Elements(t)
402
403
IF( BcMode ) THEN
404
IF( ALL( MasterEntities /= Element % BoundaryInfo % Constraint ) ) CYCLE
405
ELSE
406
IF( ALL( MasterEntities /= Element % BodyId ) ) CYCLE
407
END IF
408
409
n = Element % TYPE % NumberOfNodes
410
NodeIndexes => Element % NodeIndexes
411
IF( ElementalVar ) THEN
412
DofIndexes => Element % DGIndexes
413
ELSE
414
DofIndexes => NodeIndexes
415
END IF
416
417
DO i=1,n
418
j = DofIndexes(i)
419
IF( ASSOCIATED( Var % Perm ) ) THEN
420
k = Var % Perm(j)
421
IF( k == 0 ) CYCLE
422
ELSE
423
k = j
424
END IF
425
426
IF( VisitNodeOnlyOnce ) THEN
427
IF( VisitedNode(j) ) CYCLE
428
VisitedNode(j) = .TRUE.
429
END IF
430
431
globalnode = NodeIndexes(i)
432
433
! Only compute the parallel reduction once
434
IF( ParEnv % PEs > 1 ) THEN
435
IF( Mesh % ParallelInfo % NeighbourList(globalnode) % Neighbours(1) /= ParEnv % MyPE ) CYCLE
436
END IF
437
438
IF(NoDofs <= 1) THEN
439
x = Var % Values(k)
440
ELSE
441
x = 0.0d0
442
DO l=1,NoDofs
443
x = x + Var % Values(NoDofs*(k-1)+l) ** 2
444
END DO
445
x = SQRT(x)
446
END IF
447
448
IF(.NOT. Initialized) THEN
449
Initialized = .TRUE.
450
Maximum = x
451
Minimum = x
452
AbsMaximum = x
453
AbsMinimum = x
454
END IF
455
456
sumi = sumi + 1
457
sumx = sumx + x
458
sumxx = sumxx + x*x
459
sumabsx = sumabsx + ABS( x )
460
Maximum = MAX(x,Maximum)
461
Minimum = MIN(x,Minimum)
462
IF(ABS(x) > ABS(AbsMaximum) ) AbsMaximum = x
463
IF(ABS(x) < ABS(AbsMinimum) ) AbsMinimum = x
464
END DO
465
END DO
466
467
468
sumi = ParallelReduction(sumi)
469
IF( sumi == 0 ) THEN
470
CALL Warn('ComponentNodalReduction','No active nodes to reduced!')
471
RETURN
472
END IF
473
474
475
SELECT CASE(OperName)
476
477
CASE ('sum')
478
sumx = ParallelReduction(sumx)
479
operx = sumx
480
481
CASE ('sum abs')
482
sumx = ParallelReduction(sumabsx)
483
operx = sumabsx
484
485
CASE ('min')
486
minimum = ParallelReduction(minimum,1)
487
operx = Minimum
488
489
CASE ('max')
490
maximum = ParallelReduction(maximum,2)
491
operx = Maximum
492
493
CASE ('min abs')
494
Absminimum = ParallelReduction(AbsMinimum,1)
495
operx = AbsMinimum
496
497
CASE ('max abs')
498
Absmaximum = ParallelReduction(AbsMaximum,2)
499
operx = AbsMaximum
500
501
CASE ('range')
502
minimum = ParallelReduction(minimum,1)
503
maximum = ParallelReduction(maximum,2)
504
operx = Maximum - Minimum
505
506
CASE ('mean')
507
sumx = ParallelReduction(sumx)
508
operx = sumx / sumi
509
510
CASE ('mean abs')
511
sumx = ParallelReduction(sumabsx)
512
operx = sumabsx / sumi
513
514
CASE ('variance')
515
sumx = ParallelReduction(sumx)
516
sumxx = ParallelReduction(sumxx)
517
Operx = SQRT( sumxx/sumi-(sumx*sumx)/(sumi*sumi) )
518
519
CASE DEFAULT
520
CALL Warn('ComponentNodalReduction','Unknown statistical operator!')
521
522
END SELECT
523
524
CALL Info('ComponentNodalReduction','Reduction operator finished',Level=12)
525
526
!------------------------------------------------------------------------------
527
END FUNCTION ComponentNodalReduction
528
!------------------------------------------------------------------------------
529
530
531
!------------------------------------------------------------------------------
532
!> Perform reduction from distributed data to components.
533
!> These are similar operations as the stastistical operations in SaveScalars.
534
!------------------------------------------------------------------------------
535
FUNCTION ComponentIntegralReduction(Model, Mesh, CompParams, Var, &
536
OperName, CoeffName, GotCoeff ) RESULT ( OperX )
537
!------------------------------------------------------------------------------
538
TYPE(Model_t) :: Model
539
TYPE(Mesh_t), POINTER :: Mesh
540
TYPE(ValueList_t), POINTER :: CompParams
541
TYPE(Variable_t), POINTER :: Var
542
CHARACTER(LEN=*) :: OperName, CoeffName
543
LOGICAL :: GotCoeff
544
REAL(KIND=dp) :: OperX
545
!------------------------------------------------------------------------------
546
! Local variables
547
!------------------------------------------------------------------------------
548
TYPE(Element_t), POINTER :: Element
549
INTEGER :: t, i, j, k, n, NoDofs
550
INTEGER, POINTER :: NodeIndexes(:), DofIndexes(:)
551
REAL(KIND=dp) :: SqrtElementMetric,U,V,W,S,x,Grad(3)
552
REAL(KIND=dp) :: func, CoeffAtIp, integral, vol
553
LOGICAL :: ElementalVar, Found, Stat
554
INTEGER :: PermIndexes(Model % MaxElementNodes)
555
REAL(KIND=dp) :: Basis(Model % MaxElementNodes), dBasisdx(Model % MaxElementNodes,3)
556
REAL(KIND=dp) :: Coeff(Model % MaxElementNodes)
557
TYPE(GaussIntegrationPoints_t) :: IntegStuff
558
INTEGER, POINTER :: MasterEntities(:)
559
TYPE(Nodes_t), SAVE :: ElementNodes
560
LOGICAL, SAVE :: AllocationsDone = .FALSE.
561
INTEGER :: FirstElem, LastElem
562
LOGICAL :: BcMode
563
564
565
CALL Info('ComponentIntegralReduction','Performing reduction for component: '&
566
//TRIM(ListGetString(CompParams,'Name')),Level=10)
567
568
OperX = 0.0_dp
569
vol = 0.0_dp
570
integral = 0.0_dp
571
572
BcMode = .FALSE.
573
MasterEntities => ListGetIntegerArray( CompParams,'Master Bodies',Found )
574
IF( .NOT. Found ) THEN
575
MasterEntities => ListGetIntegerArray( CompParams,'Master Boundaries',Found )
576
BcMode = .TRUE.
577
END IF
578
579
IF(.NOT. Found ) THEN
580
CALL Warn('ComponentIntegralReduction',&
581
'> Master Bodies < or > Master Boundaries < not given')
582
RETURN
583
END IF
584
585
IF( BcMode ) THEN
586
FirstElem = Mesh % NumberOfBulkElements + 1
587
LastElem = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
588
ELSE
589
FirstElem = 1
590
LastElem = Mesh % NumberOfBulkElements
591
END IF
592
593
IF(.NOT. AllocationsDone ) THEN
594
n = Model % MaxElementNodes
595
ALLOCATE( ElementNodes % x(n), ElementNodes % y(n), ElementNodes % z(n) )
596
AllocationsDone = .TRUE.
597
END IF
598
599
ElementalVar = ( Var % TYPE == Variable_on_nodes_on_elements )
600
NoDofs = Var % Dofs
601
CoeffAtIp = 1.0_dp
602
603
DO t=FirstElem,LastElem
604
Element => Mesh % Elements(t)
605
IF( BcMode ) THEN
606
IF( ALL( MasterEntities /= Element % BoundaryInfo % Constraint ) ) CYCLE
607
ELSE
608
IF( ALL( MasterEntities /= Element % BodyId ) ) CYCLE
609
END IF
610
611
n = Element % TYPE % NumberOfNodes
612
NodeIndexes => Element % NodeIndexes
613
IF( ElementalVar ) THEN
614
DofIndexes => Element % DGIndexes
615
ELSE
616
DofIndexes => NodeIndexes
617
END IF
618
619
IF( ASSOCIATED( Var % Perm ) ) THEN
620
PermIndexes = Var % Perm( DofIndexes )
621
ELSE
622
PermIndexes = DofIndexes
623
END IF
624
625
IF ( ANY(PermIndexes == 0 ) ) CYCLE
626
627
ElementNodes % x(1:n) = Mesh % Nodes % x(NodeIndexes(1:n))
628
ElementNodes % y(1:n) = Mesh % Nodes % y(NodeIndexes(1:n))
629
ElementNodes % z(1:n) = Mesh % Nodes % z(NodeIndexes(1:n))
630
631
IF(GotCoeff) THEN
632
k = ListGetInteger( Model % Bodies( Element % BodyId ) % Values, &
633
'Material', Found )
634
Coeff(1:n) = ListGetReal( Model % Materials(k) % Values, &
635
CoeffName, n, NodeIndexes(1:n) )
636
END IF
637
638
!------------------------------------------------------------------------------
639
! Numerical integration
640
!------------------------------------------------------------------------------
641
IntegStuff = GaussPoints( Element )
642
643
DO i=1,IntegStuff % n
644
U = IntegStuff % u(i)
645
V = IntegStuff % v(i)
646
W = IntegStuff % w(i)
647
!------------------------------------------------------------------------------
648
! Basis function values & derivatives at the integration point
649
!------------------------------------------------------------------------------
650
stat = ElementInfo( Element,ElementNodes,U,V,W,SqrtElementMetric, &
651
Basis,dBasisdx )
652
!------------------------------------------------------------------------------
653
! Coordinatesystem dependent info
654
!------------------------------------------------------------------------------
655
s = SqrtElementMetric * IntegStuff % s(i)
656
IF ( CurrentCoordinateSystem() /= Cartesian ) THEN
657
x = 2 * PI * SUM( ElementNodes % x(1:n)*Basis(1:n) )
658
END IF
659
660
IF( GotCoeff ) THEN
661
CoeffAtIp = SUM( Coeff(1:n) * Basis(1:n) )
662
END IF
663
vol = vol + S
664
665
SELECT CASE(OperName)
666
667
CASE ('volume')
668
integral = integral + coeffAtIp * S
669
670
CASE ('int','int mean')
671
func = SUM( Var % Values(PermIndexes(1:n)) * Basis(1:n) )
672
integral = integral + S * coeffAtIp * func
673
674
CASE ('int abs','int abs mean')
675
func = ABS( SUM( Var % Values(PermIndexes(1:n)) * Basis(1:n) ) )
676
integral = integral + S * coeffAtIp * func
677
678
CASE ('diffusive energy')
679
DO j = 1, 3
680
Grad(j) = SUM( dBasisdx(1:n,j) * Var % Values(PermIndexes(1:n) ) )
681
END DO
682
integral = integral + s * CoeffAtIp * SUM( Grad * Grad )
683
684
CASE ('convective energy')
685
func = SUM( Var % Values(PermIndexes(1:n)) * Basis(1:n) )
686
687
IF(NoDofs == 1) THEN
688
func = SUM( Var % Values(PermIndexes(1:n)) * Basis(1:n) )
689
integral = integral + s * coeffAtIp * func**2
690
ELSE
691
func = 0.0d0
692
DO j=1,NoDofs
693
func = SUM( Var % Values(NoDofs*(PermIndexes(1:n)-1)+j) * Basis(1:n) )
694
integral = integral + s * coeffAtIp * func**2
695
END DO
696
END IF
697
698
CASE ('potential energy')
699
func = SUM( Var % Values(PermIndexes(1:n)) * Basis(1:n) )
700
integral = integral + s * coeffAtIp * func
701
702
CASE DEFAULT
703
CALL Warn('ComponentIntegralReduction','Unknown operator')
704
705
END SELECT
706
707
END DO
708
709
END DO
710
711
integral = ParallelReduction( integral )
712
713
SELECT CASE(OperName)
714
715
CASE ('volume')
716
operx = integral
717
718
CASE ('int')
719
operx = integral
720
721
CASE ('int abs')
722
operx = integral
723
724
CASE ('int mean')
725
vol = ParallelReduction( vol )
726
operx = integral / vol
727
728
CASE ('int abs mean')
729
vol = ParallelReduction( vol )
730
operx = integral / vol
731
732
CASE ('diffusive energy')
733
operx = 0.5d0 * integral
734
735
CASE ('convective energy')
736
operx = 0.5d0 * integral
737
738
CASE ('potential energy')
739
operx = integral
740
741
END SELECT
742
743
CALL Info('ComponentIntegralReduction','Reduction operator finished',Level=12)
744
745
!------------------------------------------------------------------------------
746
END FUNCTION ComponentIntegralReduction
747
!------------------------------------------------------------------------------
748
749
750
!------------------------------------------------------------------------------
751
!> Each solver may include a list of dependent components that are updated
752
!> after the solver (or the nonlinear iteration related to it) has been executed.
753
!------------------------------------------------------------------------------
754
SUBROUTINE UpdateDependentComponents( ComponentList )
755
INTEGER, POINTER :: ComponentList(:)
756
757
INTEGER :: i,j,NoVar
758
CHARACTER(:), ALLOCATABLE :: OperName, VarName, CoeffName, TmpOper
759
LOGICAL :: GotVar, GotOper, GotCoeff, VectorResult
760
TYPE(ValueList_t), POINTER :: CompParams
761
REAL(KIND=dp) :: ScalarVal, VectorVal(3), Power, Voltage
762
TYPE(Variable_t), POINTER :: Var
763
764
CALL Info('UpdateDepedentComponents','Updating Components to reflect new solution',Level=6)
765
766
DO j=1,CurrentModel % NumberOfComponents
767
768
IF( ALL( ComponentList /= j ) ) CYCLE
769
770
CALL Info('UpdateDepedentComponents','Updating component: '//I2S(j))
771
CompParams => CurrentModel % Components(j) % Values
772
773
NoVar = 0
774
DO WHILE( .TRUE. )
775
NoVar = NoVar + 1
776
OperName = ListGetString( CompParams,'Operator '//I2S(NoVar), GotOper)
777
VarName = ListGetString( CompParams,'Variable '//I2S(NoVar), GotVar)
778
CoeffName = ListGetString( CompParams,'Coefficient '//I2S(NoVar), GotCoeff)
779
780
IF(.NOT. GotVar .AND. GotOper .AND. OperName == 'electric resistance') THEN
781
VarName = 'Potential'
782
GotVar = .TRUE.
783
CALL Info('UpdateDependentComponents',&
784
'Defaulting field to > Potential < for operator: '//TRIM(OperName),Level=8)
785
END IF
786
787
IF(.NOT. (GotVar .AND. GotOper ) ) EXIT
788
789
Var => VariableGet( CurrentModel % Mesh % Variables, VarName )
790
IF( .NOT. ASSOCIATED( Var ) ) THEN
791
CALL Info('UpdateDependentComponents','Variable not available: '//TRIM(VarName))
792
CYCLE
793
END IF
794
VectorResult = .FALSE.
795
796
SELECT CASE( OperName )
797
798
CASE('electric resistance')
799
IF(.NOT. GotCoeff ) THEN
800
CoeffName = 'electric conductivity'
801
GotCoeff = .TRUE.
802
END IF
803
TmpOper = 'diffusive energy'
804
Power = ComponentIntegralReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
805
TmpOper, CoeffName, GotCoeff )
806
TmpOper = 'range'
807
Voltage = ComponentNodalReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
808
TmpOper )
809
ScalarVal = Voltage**2 / Power
810
CALL ListAddConstReal( CompParams,'res: '//TRIM(OperName),ScalarVal )
811
812
CASE ('sum','sum abs','min','max','min abs','max abs','range','mean','mean abs','variance')
813
ScalarVal = ComponentNodalReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
814
OperName )
815
CALL ListAddConstReal( CompParams,'res: '//TRIM(OperName)//' '//TRIM(VarName),ScalarVal )
816
817
CASE ('volume','int','int abs','int mean','int abs mean','diffusive energy',&
818
'convective energy','potential energy')
819
ScalarVal = ComponentIntegralReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
820
OperName, CoeffName, GotCoeff )
821
CALL ListAddConstReal( CompParams,'res: '//TRIM(OperName)//' '//TRIM(VarName),ScalarVal )
822
823
CASE('force')
824
CALL ComponentNodalForceReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
825
Force = VectorVal )
826
VectorResult = .TRUE.
827
828
CASE('moment')
829
CALL ComponentNodalForceReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
830
Moment = VectorVal )
831
VectorResult = .TRUE.
832
833
CASE('torque')
834
CALL ComponentNodalForceReduction(CurrentModel, CurrentModel % Mesh, CompParams, Var, &
835
Torque = ScalarVal )
836
837
CASE DEFAULT
838
CALL Fatal('UpdateDependentComponents','Uknown operator: '//TRIM(OperName))
839
840
END SELECT
841
842
IF( VectorResult ) THEN
843
DO i=1,3
844
WRITE( Message,'(A,ES15.6)') TRIM(OperName)//': '//TRIM(VarName)//' '&
845
//I2S(i)//': ',ScalarVal
846
CALL Info('UpdateDependentComponents',Message,Level=5)
847
CALL ListAddConstReal( CompParams,'res: '//TRIM(OperName)//': '&
848
//TRIM(VarName)//' '//I2S(i),VectorVal(i) )
849
END DO
850
ELSE
851
WRITE( Message,'(A,ES15.6)') &
852
'comp '//I2S(j)//': '//TRIM(OperName)//': '//TRIM(VarName)//': ',ScalarVal
853
CALL Info('UpdateDependentComponents',Message,Level=5)
854
CALL ListAddConstReal( CurrentModel % Simulation, &
855
'res: comp '//I2S(j)//': '//TRIM(OperName)//' '//TRIM(VarName),ScalarVal )
856
END IF
857
858
END DO
859
END DO
860
861
!------------------------------------------------------------------------------
862
END SUBROUTINE UpdateDependentComponents
863
!------------------------------------------------------------------------------
864
865
!------------------------------------------------------------------------------
866
!> Given a vector field compute line integral of Stokes theorem using
867
!> some geometric heuristics.
868
!------------------------------------------------------------------------------
869
FUNCTION ComponentStokesTheorem(Model, Mesh, Vlist, AVar, Surf ) RESULT ( FL )
870
!------------------------------------------------------------------------------
871
TYPE(Model_t) :: Model
872
TYPE(Mesh_t), POINTER :: Mesh
873
TYPE(ValueList_t), POINTER :: Vlist
874
TYPE(Variable_t), POINTER :: aVar
875
LOGICAL :: Surf
876
REAL(KIND=dp) :: FL
877
878
TYPE(Matrix_t), POINTER :: NodeGraph
879
REAL(KIND=dp), POINTER :: HelperArray(:,:)
880
REAL(KIND=dp) :: Center(3), Coord(3), Coord0(3), Coord1(3), Coord2(3), &
881
Normal(3), Tangent1(3), Tangent2(3), x1, y1, x2, y2, phi, phi1, phi2, &
882
phisum, g, ssum, h, hmin, hmax, htol, hgoal, hrel, r, rmin, rmax, rtol
883
COMPLEX :: gradv(3)
884
REAL(KIND=dp) :: EdgeVector(3), ds, r2, r2min, dsmax, dphi, dphimax, &
885
ReCirc, ImCirc
886
INTEGER :: nsteps, sgn, t, i, j, k, i1, i2, j1, j2, l, &
887
lp, n0, r2ind, imax, kmax, n, InsideBC
888
INTEGER, POINTER :: NodeIndexes(:)
889
INTEGER, POINTER :: TargetBodies(:)
890
LOGICAL :: Found, SaveLoop, Debug, GotHtol, GotRtol, BCMode
891
TYPE(Element_t), POINTER :: Edge, Element
892
LOGICAL, ALLOCATABLE :: NodeActive(:),Inside(:),Outside(:)
893
CHARACTER(*), PARAMETER :: Caller = 'ComponentStokesTheorem'
894
895
TargetBodies => ListGetIntegerArray( VList,'Master Bodies',Found )
896
IF( .NOT. Found ) TargetBodies => ListGetIntegerArray( VList,'Body',Found )
897
IF( .NOT. Found ) CALL Fatal(Caller,'Stokes theorem requires > Master Bodies <')
898
899
HelperArray => ListGetConstRealArray( Vlist, 'Coil Center', UnfoundFatal = .TRUE. )
900
Center(1:3) = HelperArray(1:3,1)
901
902
HelperArray => ListGetConstRealArray( Vlist, 'Coil normal', UnfoundFatal = .TRUE. )
903
Normal(1:3) = HelperArray(1:3,1)
904
CALL TangentDirections(Normal, Tangent1, Tangent2)
905
906
n = Mesh % NumberOfNodes
907
ALLOCATE(NodeActive(n))
908
909
InsideBC = ListGetInteger( Vlist,'Inside Boundary', BCMode )
910
911
CALL SetActiveNodeSet()
912
913
IF(Surf) THEN
914
CALL Info(Caller,'Calculating cylinder integral for: '//TRIM(avar % name))
915
CALL ComputeCylinderIntegral()
916
ELSE
917
CALL Info(Caller,'Calculating line integral for: '//TRIM(avar % name))
918
CALL ComputeLineIntegral()
919
END IF
920
921
FL = ReCirc
922
923
924
CONTAINS
925
926
! Set active nodes such that when we try to create a closed circle we will only check the marked
927
! nodes. We need a quick look-up table since otherwise there is no quick way to determine whether
928
! a node is part of a suitable edge. Also compute the height of the domain for other purposes.
929
!------------------------------------------------------------------------------------------------
930
931
SUBROUTINE SetActiveNodeSet()
932
933
IF( BCMode ) THEN
934
DO t=Mesh % NumberOfBulkElements+1, &
935
Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
936
Element => Mesh % Elements(t)
937
n = Element % TYPE % NumberOfNodes
938
NodeIndexes => Element % NodeIndexes
939
940
IF(.NOT. ASSOCIATED(Element % BoundaryInfo)) CYCLE
941
IF(Element % BoundaryInfo % Constraint == 0) CYCLE
942
DO i=1,CurrentModel % NumberOfBCs
943
IF ( Element % BoundaryInfo % Constraint == CurrentModel % BCs(i) % Tag ) EXIT
944
END DO
945
IF(i /= InsideBC) CYCLE
946
NodeActive(NodeIndexes) = .TRUE.
947
END DO
948
ELSE
949
! No BC elements given, initialize the list of candidate nodes using
950
! an intersection between inside and outside nodes.
951
ALLOCATE(Inside(n), Outside(n))
952
Inside = .FALSE.
953
Outside = .FALSE.
954
955
DO t=1,Mesh % NumberOfBulkElements
956
Element => Mesh % Elements(t)
957
IF( ANY(TargetBodies == Element % BodyId) ) THEN
958
Inside(Element % NodeIndexes) = .TRUE.
959
ELSE
960
Outside(Element % NodeIndexes) = .TRUE.
961
END IF
962
END DO
963
964
! We make the stokes theorem on the interface
965
NodeActive = Inside .AND. Outside
966
DEALLOCATE(Inside,Outside)
967
END IF
968
969
n = COUNT(NodeActive)
970
CALL Info(Caller,'Active nodes for edge patch candidates: '//I2S(n))
971
972
! We may limit the active set of nodes through which path may be drawn.
973
! The idea could be to choose only the upper or lower surface of the coil.
974
htol = ListGetConstReal( Vlist,'Flux linkage height tolerance',GotHtol )
975
rtol = ListGetConstReal( Vlist,'Flux linkage radius tolerance',GotRtol )
976
977
! We always make one round to get the bounding box!
978
IF(.TRUE.) THEN
979
hgoal = ListGetConstReal( Vlist,'Flux linkage relative height',Found )
980
hmin = HUGE(hmin)
981
hmax = -HUGE(hmax)
982
rmin = HUGE(rmin)
983
rmax = -HUGE(rmax)
984
DO k=1,2
985
DO i=1,Mesh % NumberOfNodes
986
IF(.NOT. NodeActive(i)) CYCLE
987
Coord(1) = Mesh % Nodes % x(i)
988
Coord(2) = Mesh % Nodes % y(i)
989
Coord(3) = Mesh % Nodes % z(i)
990
h = SUM(Coord*Normal)
991
r = SQRT(SUM(Coord*Tangent1)**2 + SUM(Coord*Tangent2)**2)
992
993
IF(k==1) THEN
994
hmin = MIN(h,hmin)
995
hmax = MAX(h,hmax)
996
rmin = MIN(r,rmin)
997
rmax = MAX(r,rmax)
998
ELSE
999
IF(GotHtol) THEN
1000
hrel = (h-hmin)/(hmax-hmin)
1001
IF( ABS(hrel-hgoal) > htol ) NodeActive(i) = .FALSE.
1002
END IF
1003
IF(GotRtol) THEN
1004
IF( ABS(r-rmin)/rmin > rtol) NodeActive(i) = .FALSE.
1005
END IF
1006
END IF
1007
END DO
1008
1009
IF(k==1) THEN
1010
IF( ParEnv % PEs > 1 ) THEN
1011
hmin = ParallelReduction(hmin,1)
1012
hmax = ParallelReduction(hmax,2)
1013
rmin = ParallelReduction(rmin,1)
1014
rmax = ParallelReduction(rmax,2)
1015
END IF
1016
IF(InfoActive(20)) THEN
1017
PRINT *,'Height interval in n-t system: ',hmin,hmax
1018
PRINT *,'Radius interval in n-t system: ',hmin,hmax
1019
END IF
1020
END IF
1021
1022
IF(.NOT. (GotHTol .OR. GotRTol ) ) EXIT
1023
END DO
1024
END IF
1025
1026
IF(GotHTol .OR. GotRTol ) THEN
1027
n = COUNT(NodeActive)
1028
CALL Info(Caller,'Active nodes after tolerance check (out of '&
1029
//I2S(SIZE(NodeActive))//'): '//I2S(n))
1030
END IF
1031
1032
END SUBROUTINE SetActiveNodeSet
1033
1034
1035
! This compute a line integral. No line may exist so we go through a line created on-the-fly from
1036
! existing node-to-node connections. We rotate through the axis until 360 degrees are passed
1037
! going always to the node among the candidate nodes that is as efficient as possible in terms of angle.
1038
! Unfortunately this logic may fail if the circle is not a true circle.
1039
!------------------------------------------------------------------------------------------------------
1040
SUBROUTINE ComputeLineIntegral()
1041
1042
INTEGER :: WhoActive, PrevWhoActive, NoStat, kprev(3)
1043
LOGICAL :: FileOpen
1044
REAL(KIND=dp) :: r2min_par
1045
CHARACTER(:), ALLOCATABLE :: str
1046
1047
FileOpen = .FALSE.
1048
NoStat = 0
1049
1050
1051
! Create a graph for node-to-edge connectivity
1052
!----------------------------------------------
1053
NodeGraph => AllocateMatrix()
1054
NodeGraph % FORMAT = MATRIX_LIST
1055
DO i = Mesh % NumberOfEdges, 1, -1
1056
Edge => Mesh % Edges(i)
1057
IF(ALL(NodeActive(Edge % NodeIndexes))) THEN
1058
DO j=1, Edge % TYPE % NumberOfNodes
1059
CALL List_AddToMatrixElement( NodeGraph % ListMatrix,Edge % NodeIndexes(j),i,1.0_dp )
1060
END DO
1061
END IF
1062
END DO
1063
CALL List_ToCRSMatrix(NodeGraph)
1064
WRITE(Message,'(A,ES12.3)') 'Nonzeros per row NodeGraph:',&
1065
1.0_dp * SIZE(NodeGraph % Values) / NodeGraph % NumberOfRows
1066
CALL Info(Caller, Message, Level=10)
1067
1068
n0 = Mesh % NumberOfNodes
1069
1070
! Find the node closest to the origin
1071
r2min = HUGE(r2min)
1072
r2ind = 0
1073
DO i=1,Mesh % NumberOfNodes
1074
IF(.NOT. NodeActive(i)) CYCLE
1075
Coord(1) = Mesh % Nodes % x(i)
1076
Coord(2) = Mesh % Nodes % y(i)
1077
Coord(3) = Mesh % Nodes % z(i)
1078
r2 = SUM((Coord-Center)**2)
1079
IF(r2 < r2min) THEN
1080
r2min = r2
1081
r2ind = i
1082
Coord0 = Coord
1083
END IF
1084
END DO
1085
1086
WRITE(Message,'(A,ES10.3)') 'Minimum distance (node '//I2S(r2ind)//'):',SQRT(r2min)
1087
CALL Info(Caller,Message,Level=7)
1088
1089
Debug = .FALSE.
1090
IF( Debug ) THEN
1091
PRINT *,'Center:',Center
1092
PRINT *,'Normal:',Normal
1093
PRINT *,'Tangent1:',Tangent1
1094
PRINT *,'Tangent2:',Tangent2
1095
PRINT *,'Coord0:',Coord0
1096
END IF
1097
1098
phisum = 0.0_dp
1099
nsteps = 0
1100
ReCirc = 0.0_dp
1101
ImCirc = 0.0_dp
1102
ssum = 0.0_dp
1103
1104
Coord1 = Coord0
1105
i1 = r2ind
1106
1107
! We have to assume that one partition finds the entire circle.
1108
! If this is done in many pieces we should generate an algo that passed on the
1109
! process on a joint node that turns out to be otherwise the end of the path.
1110
IF(ParEnv % PEs > 1) THEN
1111
WhoActive = -1
1112
r2min_par = ParallelReduction(r2min,1)
1113
IF( ABS(r2min_par - r2min) < TINY(r2min)) THEN
1114
WhoActive = ParEnv % MyPe
1115
END IF
1116
1117
! We start from this partition.
1118
! I.e. the largest partition index with the minimum distance.
1119
WhoActive = ParallelReduction(WhoActive,2)
1120
1121
! We need Coord0 always in parallel in each partition.
1122
IF(WhoActive /= ParEnv % MyPe ) THEN
1123
Coord0 = -HUGE(Coord0)
1124
r2ind = 0
1125
END IF
1126
DO k=1,3
1127
Coord0(k) = ParallelReduction(Coord0(k),2)
1128
END DO
1129
ELSE
1130
WhoActive = ParEnv % Mype
1131
END IF
1132
PrevWhoActive = WhoActive
1133
1134
10 IF( WhoActive == ParEnv % MyPe ) THEN
1135
x1 = SUM(Coord1*Tangent1)
1136
y1 = SUM(Coord1*Tangent2)
1137
phi1 = (180.0_dp/PI) * ATAN2(y1,x1)
1138
str = ListGetString( Vlist,'Line Integral File',SaveLoop )
1139
ELSE
1140
Coord1 = -HUGE(Coord1)
1141
SaveLoop = .FALSE.
1142
END IF
1143
1144
IF( ParEnv % PEs > 1) THEN
1145
! If this is not active partition go to wait for the active one for further instructions.
1146
IF( WhoActive /= ParEnv % MyPe ) GOTO 20
1147
END IF
1148
1149
IF( SaveLoop ) THEN
1150
! Only open the file once!
1151
IF(.NOT. FileOpen ) THEN
1152
IF( ParEnv % PEs > 1 ) THEN
1153
OPEN (10, FILE=TRIM(str)//'_'//I2S(ParEnv % MyPe) )
1154
ELSE
1155
OPEN (10, FILE=str )
1156
END IF
1157
FileOpen = .TRUE.
1158
END IF
1159
WRITE(10,*) nsteps, Coord1, phi1, phisum, ssum, ReCirc
1160
END IF
1161
1162
kprev = 0
1163
DO WHILE(.TRUE.)
1164
dsmax = -EPSILON(dsmax)
1165
kprev(2:3) = kprev(1:2)
1166
kprev(1) = kmax
1167
kmax = 0
1168
1169
! Among the edges related to node "i1" find the one that has makes us further in
1170
! minimizing the distance.
1171
DO j = NodeGraph % Rows(i1),NodeGraph % Rows(i1+1)-1
1172
k = NodeGraph % Cols(j)
1173
1174
! Do not use any of the previous node again!
1175
IF(ANY(k==kprev)) CYCLE
1176
Edge => Mesh % Edges(k)
1177
NodeIndexes => Edge % NodeIndexes
1178
1179
IF(NodeIndexes(1) == i1) THEN
1180
i2 = NodeIndexes(2)
1181
ELSE IF(NodeIndexes(2) == i1) THEN
1182
i2 = NodeIndexes(1)
1183
ELSE
1184
CALL Fatal(Caller,'Either node index in edge should be i1!')
1185
END IF
1186
IF(.NOT. NodeActive(i2)) CYCLE
1187
1188
Coord(1) = Mesh % Nodes % x(i2)
1189
Coord(2) = Mesh % Nodes % y(i2)
1190
Coord(3) = Mesh % Nodes % z(i2)
1191
1192
x2 = SUM(Coord*Tangent1)
1193
y2 = SUM(Coord*Tangent2)
1194
phi = (180.0_dp/PI) * ATAN2(y2,x2)
1195
1196
dphi = phi1-phi
1197
1198
! Deal with the discontinuity of angle around +/-180 deg
1199
IF(dphi < -180.0_dp) dphi = dphi+360.0_dp
1200
IF(dphi > 180.0_dp) dphi = dphi-360.0_dp
1201
1202
IF(phisum < 315.0_dp ) THEN
1203
! This measure takes the shortest route at least in some case more robustly than some others.
1204
ds = dphi / SUM((Coord1-Coord)**2)
1205
ELSE
1206
! After we only have ~45 degs left we find the node which approaches the starting node
1207
! as well as possible.
1208
ds = SUM( (Coord1-Coord0)**2 - (Coord-Coord0)**2 ) / SUM((Coord1-Coord)**2)
1209
END IF
1210
1211
IF( Debug ) THEN
1212
PRINT *,'kmax:',ds > dsmax, k,i2,ds,dphi,phi,Coord
1213
END IF
1214
1215
! Have we found a better edge candidate? Memorize that!
1216
IF( ds >= dsmax ) THEN
1217
kmax = k
1218
imax = i2
1219
dsmax = ds
1220
Coord2 = Coord
1221
dphimax = dphi
1222
phi2 = phi
1223
END IF
1224
END DO
1225
1226
! When no way to get further we are done!
1227
IF( kmax == 0 ) THEN
1228
IF( ParEnv % PEs > 1 ) THEN
1229
BLOCK
1230
INTEGER, POINTER :: Neig(:)
1231
Neig => Mesh % ParallelInfo % NeighbourList(i1) % Neighbours
1232
IF( SIZE(Neig) > 1 ) THEN
1233
!PRINT *,'Swapping to partition:',Neig
1234
IF( WhoActive == Neig(1) ) THEN
1235
WhoActive = Neig(2)
1236
ELSE
1237
WhoActive = Neig(1)
1238
END IF
1239
EXIT
1240
END IF
1241
END BLOCK
1242
END IF
1243
1244
!PRINT *,'Cands:',i1,NodeGraph % Rows(i1),NodeGraph % Rows(i1+1)-1
1245
NoStat = 1
1246
CALL Warn(Caller,'We had to stop because no route was found!')
1247
EXIT
1248
END IF
1249
1250
nsteps = nsteps + 1
1251
1252
! Edge that goes to the minimum value
1253
k = kmax
1254
Edge => Mesh % Edges(k)
1255
i2 = imax
1256
1257
EdgeVector = Coord2 - Coord1
1258
1259
! Integration length to check optimality of the path
1260
ds = SQRT(SUM(EdgeVector**2))
1261
ssum = ssum + ds
1262
1263
IF( MODULO(avar % dofs,3) == 0 ) THEN
1264
! Integral over nodal field using the mean value
1265
j1 = avar % Perm(i1)
1266
j2 = avar % Perm(i2)
1267
IF(j1==0 .OR. j2==0) CALL Fatal(Caller,'Nodal field missing on path!')
1268
DO k=1,3
1269
IF( avar % dofs == 3 ) THEN
1270
gradv(k) = ( avar % Values(3*(j1-1)+k) + avar % Values(3*(j2-1)+k) ) / 2
1271
ELSE
1272
gradv(k) = CMPLX(avar % Values(6*(j1-1)+k) + avar % Values(6*(j2-1)+k),&
1273
avar % Values(6*(j1-1)+3+k) + avar % Values(6*(j2-1)+3+k), KIND=dp ) / 2
1274
END IF
1275
END DO
1276
ReCirc = ReCirc + REAL(SUM(gradv*EdgeVector))
1277
ImCirc = ImCirc + AIMAG(SUM(gradv*EdgeVector))
1278
ELSE
1279
! Integral over edge field.
1280
1281
! Check the sign if the direction based on global edge direction rules
1282
! If we do the path integral in the wrong direction compared to definition of edge switch the sign
1283
sgn = 1
1284
IF( ParEnv % PEs > 1 ) THEN
1285
i1 = Mesh % ParallelInfo % GlobalDOFs(i1)
1286
i2 = Mesh % ParallelInfo % GlobalDOFs(i2)
1287
END IF
1288
!IF(XOR(Edge % NodeIndexes(1) /= i1, i1 < i2) ) sgn = -sgn
1289
IF(i1 > i2) sgn = -1
1290
1291
j = avar % Perm(n0 + k)
1292
IF( j==0) CALL Fatal(Caller,'Edge field missing on path!')
1293
1294
IF( avar % dofs == 1 ) THEN
1295
ReCirc = ReCirc + sgn * avar % Values(j)
1296
ELSE
1297
ReCirc = ReCirc + sgn * avar % Values(2*j-1)
1298
ImCirc = ImCirc + sgn * avar % Values(2*j)
1299
END IF
1300
END IF
1301
1302
! Do not visit this point a 2nd time!
1303
!NodeActive(i1) = .FALSE.
1304
1305
! Continue from the end point
1306
i1 = imax
1307
Coord1 = Coord2
1308
phi1 = phi2
1309
1310
phisum = phisum + dphimax
1311
1312
IF(Debug) PRINT *,'phisum:',nsteps,Coord1,phisum,ssum,phi1,ReCirc
1313
1314
IF( ABS(phisum) > 720.0_dp ) THEN
1315
CALL Fatal(Caller,'We circled twice around!?')
1316
END IF
1317
1318
IF(SaveLoop) WRITE(10,*) nsteps, Coord1, phi1, phisum, ssum, ReCirc
1319
1320
! We have come home to roost!
1321
IF(i1 == r2ind) THEN
1322
CALL Info(Caller,'Integration over full loop finished!')
1323
EXIT
1324
END IF
1325
END DO
1326
1327
20 IF( ParEnv % PEs > 1 ) THEN
1328
phisum = ParallelReduction(phisum)
1329
!PRINT *,'phisum:',phisum, ParEnv % MyPe
1330
!PRINT *,'whoactive:',whoactive, prevwhoactive, ParEnv % Mype
1331
1332
! Find out the new active partition. It has been set only in the previous active partition.
1333
IF(PrevWhoActive /= ParEnv % MyPe) THEN
1334
WhoActive = -1
1335
END IF
1336
WhoActive = ParallelReduction(WhoActive,2)
1337
1338
!PRINT *,'New whoactive:',WhoActive, PrevWhoActive, ParEnv % MyPe
1339
1340
IF(WhoActive /= PrevWhoActive) THEN
1341
PrevWhoActive = WhoActive
1342
1343
! Find the new starting node.
1344
DO k=1,3
1345
Coord1(k) = ParallelReduction(Coord1(k),2)
1346
END DO
1347
1348
! For the new active partition find a new starting point.
1349
IF( WhoActive == ParEnv % MyPe ) THEN
1350
r2min = HUGE(r2min)
1351
i1 = 0
1352
DO i=1,Mesh % NumberOfNodes
1353
IF(.NOT. NodeActive(i)) CYCLE
1354
Coord(1) = Mesh % Nodes % x(i)
1355
Coord(2) = Mesh % Nodes % y(i)
1356
Coord(3) = Mesh % Nodes % z(i)
1357
r2 = SUM((Coord-Coord1)**2)
1358
IF(r2 < r2min) THEN
1359
r2min = r2
1360
i1 = i
1361
END IF
1362
END DO
1363
!PRINT *,'New starting point:',i1,SQRT(r2min)
1364
ELSE
1365
! Have the counter active only in the active partition.
1366
phisum = 0.0_dp
1367
END IF
1368
1369
!PRINT *,'WhoIsActice:',ParEnv % MyPe, WhoActive, WhoActive == ParEnv % MyPe
1370
GOTO 10
1371
END IF
1372
END IF
1373
1374
IF(ParEnv % PEs > 1 ) THEN
1375
ReCirc = ParallelReduction(ReCirc)
1376
ImCirc = ParallelReduction(ImCirc)
1377
ssum = ParallelReduction(ssum)
1378
END IF
1379
1380
IF(InfoActive(20)) THEN
1381
PRINT *,'PathIntegralLine:',ParEnv % MyPe, avar % dofs, targetbodies, nsteps, phisum, ssum, ReCirc
1382
END IF
1383
1384
NoStat = ParallelReduction(NoStat)
1385
IF(Nostat > 0 ) THEN
1386
CALL Info(Caller,'Setting values to zero because of issues!')
1387
ReCirc = 0.0_dp
1388
ImCirc = 0.0_dp
1389
END IF
1390
1391
CALL FreeMatrix(NodeGraph)
1392
IF(FileOpen) CLOSE(10)
1393
1394
END SUBROUTINE ComputeLineIntegral
1395
1396
1397
SUBROUTINE ComputeCylinderIntegral()
1398
1399
TYPE(Nodes_t), SAVE :: ElementNodes
1400
LOGICAL :: AllocationsDone = .FALSE.
1401
TYPE(GaussIntegrationPoints_t) :: IP
1402
REAL(KIND=dp) :: detJ, Area, s, TestVec(3), ParentNormal(3), Coeff
1403
TYPE(ValueList_t), POINTER :: Params
1404
LOGICAL :: Stat, PiolaVersion, EdgeBasis
1405
INTEGER :: np, nd, EdgeBasisDegree, tmin, tmax
1406
INTEGER, POINTER, SAVE :: Indexes(:)
1407
TYPE(Element_t), POINTER :: Element, Parent, pElem
1408
REAL(KIND=dp), POINTER, SAVE :: Basis(:), SOL(:,:), WBasis(:,:), dBasisdx(:,:), RotWBasis(:,:)
1409
1410
CALL Info(Caller,'Estimating line integral from surface integral!')
1411
1412
IF(.NOT. AllocationsDone ) THEN
1413
n = 2*Model % MaxElementNodes
1414
ALLOCATE( ElementNodes % x(n), ElementNodes % y(n), ElementNodes % z(n), &
1415
Basis(n), dBasisdx(n,3), WBasis(n,3), RotWBasis(n,3), SOL(6,n), Indexes(n) )
1416
AllocationsDone = .TRUE.
1417
END IF
1418
1419
Area = 0.0_dp
1420
ReCirc = 0.0_dp
1421
ImCirc = 0.0_dp
1422
1423
EdgeBasis = .FALSE.
1424
IF(avar % dofs <= 2) THEN
1425
EdgeBasis = .TRUE.
1426
Params => avar % Solver % Values
1427
CALL EdgeElementStyle(avar % Solver % Values, PiolaVersion, BasisDegree = EdgeBasisDegree )
1428
END IF
1429
1430
IF( BCMode ) THEN
1431
tmin = Mesh % NumberOfBulkElements + 1
1432
tmax = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
1433
ELSE
1434
tmin = 1
1435
tmax = Mesh % NumberOfFaces
1436
END IF
1437
1438
!PRINT *,'Center:',Center
1439
!PRINT *,'Normal:',Normal
1440
!PRINT *,'BCMode:',BCMode,tmin,tmax,InsideBC
1441
1442
DO t=tmin, tmax
1443
IF(BCMode) THEN
1444
Element => Mesh % Elements(t)
1445
ELSE
1446
Element => Mesh % Faces(t)
1447
END IF
1448
n = Element % TYPE % NumberOfNodes
1449
NodeIndexes => Element % NodeIndexes
1450
1451
IF( BCMode ) THEN
1452
IF(.NOT. ASSOCIATED(Element % BoundaryInfo)) CYCLE
1453
IF(Element % BoundaryInfo % Constraint == 0) CYCLE
1454
DO i=1,CurrentModel % NumberOfBCs
1455
IF ( Element % BoundaryInfo % Constraint == CurrentModel % BCs(i) % Tag ) EXIT
1456
END DO
1457
IF(i /= InsideBC) CYCLE
1458
ELSE
1459
IF(.NOT. ALL(NodeActive(NodeIndexes))) CYCLE
1460
END IF
1461
1462
! Check that we have a parent that is on the boundary.
1463
k = 0
1464
DO i=1,2
1465
IF(i==1) THEN
1466
pElem => Element % BoundaryInfo % Left
1467
ELSE
1468
pElem => Element % BoundaryInfo % Right
1469
END IF
1470
IF(ASSOCIATED(pElem)) THEN
1471
IF( ANY( TargetBodies == pElem % BodyId ) ) THEN
1472
k=k+1
1473
Parent => pElem
1474
END IF
1475
END IF
1476
END DO
1477
IF(k/=1) CYCLE
1478
1479
ElementNodes % x(1:n) = Mesh % Nodes % x(NodeIndexes(1:n))
1480
ElementNodes % y(1:n) = Mesh % Nodes % y(NodeIndexes(1:n))
1481
ElementNodes % z(1:n) = Mesh % Nodes % z(NodeIndexes(1:n))
1482
1483
! Normal pointing out of the parent
1484
ParentNormal = NormalVector(Element,ElementNodes,Parent=Parent)
1485
! Normal just in xy-plane
1486
ParentNormal = ParentNormal - SUM(Normal*ParentNormal)*Normal
1487
1488
! Center of face element
1489
Coord1(1) = SUM(ElementNodes % x(1:n)) / n
1490
Coord1(2) = SUM(ElementNodes % y(1:n)) / n
1491
Coord1(3) = SUM(ElementNodes % z(1:n)) / n
1492
! Face element direction vector in xy-plane
1493
Coord1 = Coord1 - Center
1494
Coord1 = Coord1 - SUM(Normal*Coord1)*Normal
1495
! Normalize such that |Coord1| == 1.
1496
ds = SQRT(SUM(Coord1**2))
1497
IF(ds > EPSILON(ds)) Coord1 = Coord1 / ds
1498
1499
! If we are not pointing inwards at all then skip the face element.
1500
Coeff = -SUM(ParentNormal * Coord1)
1501
IF(Coeff < EPSILON(Coeff) ) CYCLE
1502
1503
! Find the maximum distance edge on the boundary in the local coordinates.
1504
! This edge is oriented with the surface having the desired direction.
1505
dsmax = -HUGE(dsmax)
1506
DO i=1,n
1507
Coord1(1) = ElementNodes % x(i)
1508
Coord1(2) = ElementNodes % y(i)
1509
Coord1(3) = ElementNodes % z(i)
1510
Coord1 = Coord1 - Center
1511
Coord1 = Coord1 - SUM(Normal*Coord1)*Normal
1512
DO j=i+1,n
1513
Coord2(1) = ElementNodes % x(j)
1514
Coord2(2) = ElementNodes % y(j)
1515
Coord2(3) = ElementNodes % z(j)
1516
Coord2 = Coord2 - Center
1517
Coord2 = Coord2 - SUM(Normal*Coord2)*Normal
1518
ds = SUM((Coord1-Coord2)**2)
1519
IF(ds > dsmax ) THEN
1520
dsmax = ds
1521
EdgeVector = Coord1-Coord2
1522
END IF
1523
END DO
1524
END DO
1525
EdgeVector = EdgeVector / SQRT(SUM(EdgeVector**2))
1526
1527
TestVec = CrossProduct(Normal,Coord1)
1528
IF(SUM(TestVec*EdgeVector) > 0.0_dp) EdgeVector = -EdgeVector
1529
1530
IF( EdgeBasis ) THEN
1531
nd = mGetElementDofs( Indexes, Uelement = Element, USolver = avar % Solver )
1532
np = COUNT(Indexes(1:nd) <= Mesh % NumberOfNodes)
1533
IP = GaussPoints(Element, EdgeBasis=.TRUE., PReferenceElement=PiolaVersion, &
1534
EdgeBasisDegree=EdgeBasisDegree)
1535
ELSE
1536
Indexes(1:n) = NodeIndexes(1:n)
1537
nd = n
1538
IP = GaussPoints( Element )
1539
END IF
1540
1541
DO i=1,avar % dofs
1542
SOL(i,1:nd) = avar % values(avar % dofs*(avar % Perm(Indexes(1:nd))-1)+i)
1543
END DO
1544
1545
!------------------------------------------------------------------------------
1546
! Numerical integration
1547
!------------------------------------------------------------------------------
1548
1549
DO l=1,IP % n
1550
IF(.NOT. EdgeBasis) THEN
1551
stat = ElementInfo( Element,ElementNodes,&
1552
IP % U(l),IP % V(l),IP % W(l), DetJ, Basis )
1553
ELSE
1554
stat = ElementInfo( Element, ElementNodes, IP % U(l), IP % V(l), &
1555
IP % W(l), detJ, Basis, dBasisdx, EdgeBasis = WBasis, &
1556
RotBasis = RotWBasis, USolver = avar % Solver )
1557
END IF
1558
1559
! stat = EdgeElementInfo(Element, ElementNodes, IP % U(l), IP % V(l), IP % W(l), &
1560
! DetF = DetJ, Basis = Basis, EdgeBasis = WBasis, dBasisdx = dBasisdx, &
1561
! BasisDegree = EdgeBasisDegree, ApplyPiolaTransform = .TRUE.)
1562
! ELSE
1563
! stat = ElementInfo(Element, ElementNodes, IP % U(l), IP % V(l), IP % W(l), &
1564
! detJ, Basis, dBasisdx)
1565
! CALL GetEdgeBasis(Element, WBasis, RotWBasis, Basis, dBasisdx)
1566
! END IF
1567
1568
s = Coeff * DetJ * IP % s(l)
1569
Area = Area + S
1570
1571
SELECT CASE( avar % dofs )
1572
CASE( 1 )
1573
gradv = MATMUL(SOL(1,np+1:nd), WBasis(1:nd-np,:))
1574
CASE( 2 )
1575
gradv = CMPLX( MATMUL(SOL(1,np+1:nd), WBasis(1:nd-np,:)), &
1576
MATMUL(SOL(2,np+1:nd), WBasis(1:nd-np,:)), KIND=dp)
1577
CASE( 3 )
1578
gradv = MATMUL(SOL(1:3,1:n),Basis(1:n))
1579
CASE( 6 )
1580
DO i=1,3
1581
gradv(i) = CMPLX( SUM(SOL(2*i-1,1:n)*Basis(1:n)), &
1582
SUM(SOL(2*i,1:n)*Basis(1:n)), KIND=dp )
1583
END DO
1584
END SELECT
1585
1586
ReCirc = ReCirc + s * REAL(SUM(gradv*EdgeVector))
1587
ImCirc = ImCirc + s * AIMAG(SUM(gradv*EdgeVector))
1588
END DO
1589
END DO
1590
1591
! Sum up in parallel.
1592
Area = ParallelReduction(Area)
1593
ReCirc = ParallelReduction(ReCirc)
1594
ImCirc = ParallelReduction(ImCirc)
1595
1596
! Move from surface integral to line integral correspondent by dividing with the height.
1597
ReCirc = ReCirc / (hmax-hmin)
1598
ImCirc = ImCirc / (hmax-hmin)
1599
1600
IF(InfoActive(20)) THEN
1601
PRINT *,'PathIntegralCyl:',avar % dofs, targetbodies, area, &
1602
area/((hmax-hmin)*2*PI), ReCirc, ImCirc
1603
END IF
1604
1605
END SUBROUTINE ComputeCylinderIntegral
1606
1607
END FUNCTION ComponentStokesTheorem
1608
1609
1610
!------------------------------------------------------------------------------
1611
!> Given vector potential and current density compute the energy in the coil.
1612
!> This is actually not energy, but twice the energy, since the values are
1613
!> used to computed inductance matrix.
1614
!------------------------------------------------------------------------------
1615
FUNCTION ComponentCoilEnergy(Model, Mesh, MasterEntities, AVar, CVar, BCMode ) RESULT ( AIintRe )
1616
!------------------------------------------------------------------------------
1617
TYPE(Model_t) :: Model
1618
TYPE(Mesh_t), POINTER :: Mesh
1619
INTEGER, POINTER :: MasterEntities(:)
1620
TYPE(Variable_t), POINTER :: AVar, CVar
1621
LOGICAL, OPTIONAL :: BCMode
1622
REAL(KIND=dp) :: AIintRe
1623
!------------------------------------------------------------------------------
1624
! Local variables
1625
!------------------------------------------------------------------------------
1626
TYPE(Element_t), POINTER :: Element
1627
INTEGER :: t, i, j, k, l, n, np, nd, EdgeBasisDegree, t1, t2
1628
REAL(KIND=dp) :: volume
1629
LOGICAL :: Found
1630
LOGICAL :: Stat, PiolaVersion, EdgeBasis, DoBCs
1631
COMPLEX(KIND=dp) :: AIint
1632
CHARACTER(LEN=MAX_NAME_LEN) :: str
1633
CHARACTER(*), PARAMETER :: Caller = 'ComponentCoilEnergy'
1634
1635
IF(.NOT. ASSOCIATED(MasterEntities)) THEN
1636
CALL Fatal(Caller,'"MasterEntities" not associated!')
1637
END IF
1638
1639
DoBCs = .FALSE.
1640
IF(PRESENT(BcMode)) DoBCs = BCMode
1641
1642
str = I2S(MasterEntities(1))
1643
DO i=2,SIZE(MasterEntities)
1644
str = TRIM(str)//' '//I2S(MasterEntities(i))
1645
END DO
1646
IF( DoBCs ) THEN
1647
CALL Info(Caller,'Performing reduction for BCs: '//TRIM(str),Level=10)
1648
ELSE
1649
CALL Info(Caller,'Performing reduction for bodies: '//TRIM(str),Level=10)
1650
END IF
1651
1652
IF( DoBCs ) THEN
1653
t1 = Mesh % NumberOfBulkElements + 1
1654
t2 = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
1655
ELSE
1656
t1 = 1
1657
t2 = Mesh % NumberOfBulkElements
1658
END IF
1659
1660
EdgeBasis = .FALSE.
1661
IF(avar % dofs <= 2) THEN
1662
EdgeBasis = .TRUE.
1663
CALL EdgeElementStyle(avar % Solver % Values, PiolaVersion, BasisDegree = EdgeBasisDegree )
1664
END IF
1665
1666
AIint = 0.0_dp
1667
volume = 0.0_dp
1668
1669
IF( CVar % Dofs /= 3 ) THEN
1670
CALL Fatal(Caller,'Expecting 3 components for current density!')
1671
END IF
1672
IF( ALL([1,2,3,6] /= AVar % Dofs) ) THEN
1673
CALL Fatal(Caller,'Expecting 1,2,3 or 6 components for vector potential!')
1674
END IF
1675
1676
DO t=t1, t2
1677
Element => Mesh % Elements(t)
1678
IF( DoBCs ) THEN
1679
IF( ALL( MasterEntities /= Element % BoundaryInfo % Constraint ) ) CYCLE
1680
ELSE
1681
IF( ALL( MasterEntities /= Element % BodyId ) ) CYCLE
1682
END IF
1683
CALL LocalIntegElem()
1684
END DO
1685
1686
AIint = ParallelReduction( AIint )
1687
Volume = ParallelReduction( volume )
1688
1689
!PRINT *,'AiInit:',AiInt,Volume
1690
1691
AIIntRe = REAL(AIint)
1692
1693
CALL Info(Caller,'Reduction operator finished',Level=12)
1694
1695
CONTAINS
1696
1697
SUBROUTINE LocalIntegElem()
1698
1699
TYPE(GaussIntegrationPoints_t) :: IP
1700
TYPE(Nodes_t), SAVE :: ElementNodes
1701
LOGICAL, SAVE :: AllocationsDone = .FALSE.
1702
INTEGER, POINTER, SAVE :: EdgeIndexes(:)
1703
INTEGER, POINTER :: NodeIndexes(:), pIndexes(:)
1704
REAL(KIND=dp) :: DetJ,S,Cip(3)
1705
COMPLEX(KIND=dp) :: Aip(3)
1706
REAL(KIND=dp), POINTER, SAVE :: Basis(:), WBasis(:,:), dBasisdx(:,:), RotWBasis(:,:), &
1707
Aelem(:,:), Celem(:,:)
1708
1709
IF(.NOT. AllocationsDone ) THEN
1710
n = 2*Model % MaxElementNodes
1711
ALLOCATE( ElementNodes % x(n), ElementNodes % y(n), ElementNodes % z(n), &
1712
Basis(n), dBasisdx(n,3), WBasis(n,3), RotWBasis(n,3), Aelem(6,n), Celem(3,n), EdgeIndexes(n) )
1713
AllocationsDone = .TRUE.
1714
END IF
1715
1716
1717
n = Element % TYPE % NumberOfNodes
1718
NodeIndexes => Element % NodeIndexes
1719
1720
ElementNodes % x(1:n) = Mesh % Nodes % x(NodeIndexes(1:n))
1721
ElementNodes % y(1:n) = Mesh % Nodes % y(NodeIndexes(1:n))
1722
ElementNodes % z(1:n) = Mesh % Nodes % z(NodeIndexes(1:n))
1723
1724
IF( EdgeBasis ) THEN
1725
nd = mGetElementDofs( EdgeIndexes, Uelement = Element, USolver = avar % Solver )
1726
np = COUNT(EdgeIndexes(1:nd) <= Mesh % NumberOfNodes)
1727
IP = GaussPoints(Element, EdgeBasis=.TRUE., PReferenceElement=PiolaVersion, &
1728
EdgeBasisDegree=EdgeBasisDegree)
1729
pIndexes => EdgeIndexes
1730
ELSE
1731
nd = n
1732
IP = GaussPoints( Element )
1733
pIndexes => NodeIndexes
1734
END IF
1735
1736
! Edge or nodal values of vector potential
1737
DO i=1,avar % dofs
1738
Aelem(i,1:nd) = avar % values(avar % dofs*(avar % Perm(pIndexes(1:nd))-1)+i)
1739
END DO
1740
1741
! Nodal values of current density
1742
IF( cvar % TYPE == Variable_on_nodes_on_elements ) THEN
1743
pIndexes => Element % DGIndexes
1744
ELSE
1745
pIndexes => NodeIndexes
1746
END IF
1747
DO i=1,cvar % dofs
1748
Celem(i,1:n) = cvar % values(cvar % dofs*(cvar % Perm(pIndexes(1:n))-1)+i)
1749
END DO
1750
1751
DO l=1,IP % n
1752
IF(.NOT. EdgeBasis) THEN
1753
stat = ElementInfo( Element,ElementNodes,&
1754
IP % U(l),IP % V(l),IP % W(l), DetJ, Basis )
1755
ELSE
1756
stat = ElementInfo( Element, ElementNodes, IP % U(l), IP % V(l), &
1757
IP % W(l), detJ, Basis, dBasisdx, EdgeBasis = WBasis, &
1758
RotBasis = RotWBasis, USolver = avar % Solver )
1759
END IF
1760
1761
s = DetJ * IP % s(l)
1762
1763
! Vector potential at IP
1764
SELECT CASE( avar % dofs )
1765
CASE( 1 )
1766
Aip = MATMUL(Aelem(1,np+1:nd), WBasis(1:nd-np,:))
1767
CASE( 2 )
1768
Aip = CMPLX( MATMUL(Aelem(1,np+1:nd), WBasis(1:nd-np,:)), &
1769
MATMUL(Aelem(2,np+1:nd), WBasis(1:nd-np,:)), KIND=dp)
1770
CASE( 3 )
1771
Aip = MATMUL(Aelem(1:3,1:n),Basis(1:n))
1772
CASE( 6 )
1773
Aip = CMPLX( MATMUL(Aelem(1:5:2,1:n),Basis(1:n)), &
1774
MATMUL(Aelem(2:6:2,1:n),Basis(1:n)), KIND=dp )
1775
END SELECT
1776
1777
! Current density at IP
1778
Cip = MATMUL(Celem(1:3,1:n),Basis(1:n))
1779
1780
!PRINT *,'Ai:',s,SQRT(SUM(REAL(Aip)**2)),SQRT(SUM(AIMAG(Aip)**2)),'c',SQRT(SUM(Cip*Cip))
1781
1782
AIint = AIint + s * SUM(Aip*Cip)
1783
Volume = Volume + s
1784
END DO
1785
1786
END SUBROUTINE LocalIntegElem
1787
1788
!------------------------------------------------------------------------------
1789
END FUNCTION ComponentCoilEnergy
1790
!------------------------------------------------------------------------------
1791
1792
1793
!------------------------------------------------------------------------------
1794
!> Compute integrals for determining the S parameters.
1795
!------------------------------------------------------------------------------
1796
FUNCTION BoundaryWaveFlux(Model, Mesh, MasterEntities, Avar, InFlux, PortImp, PortBC ) &
1797
RESULT ( OutFlux )
1798
!------------------------------------------------------------------------------
1799
TYPE(Model_t) :: Model
1800
TYPE(Mesh_t), POINTER :: Mesh
1801
INTEGER, POINTER :: MasterEntities(:)
1802
TYPE(Variable_t), POINTER :: Avar
1803
COMPLEX(KIND=dp) :: OutFlux
1804
COMPLEX(KIND=dp) :: InFlux
1805
COMPLEX(KIND=dp) :: PortImp
1806
LOGICAL :: PortBC
1807
!------------------------------------------------------------------------------
1808
! Local variables
1809
!------------------------------------------------------------------------------
1810
TYPE(Element_t), POINTER :: Element
1811
INTEGER :: t, i, j, k, l, n, np, nd, EdgeBasisDegree, t1, t2, bc_id
1812
LOGICAL :: Found, InitHandles
1813
LOGICAL :: Stat, PiolaVersion, EdgeBasis, UseGaussLaw
1814
TYPE(ValueList_t), POINTER :: BC
1815
CHARACTER(LEN=MAX_NAME_LEN) :: str
1816
REAL(KIND=dp) :: area, omega
1817
COMPLEX(KIND=dp) :: int_norm, int_el, vol, curr, port_curr, trans, Zimp
1818
CHARACTER(*), PARAMETER :: Caller = 'BoundaryWaveFlux'
1819
1820
area = 0.0_dp
1821
1822
int_norm = 0.0_dp
1823
int_el = 0.0_dp
1824
1825
vol = 0.0_dp
1826
curr = 0.0_dp
1827
port_curr = 0.0_dp
1828
trans = 0.0_dp
1829
1830
IF(.NOT. ASSOCIATED(MasterEntities)) THEN
1831
CALL Fatal(Caller,'"MasterEntities" not associated!')
1832
END IF
1833
1834
str = I2S(MasterEntities(1))
1835
DO i=2,SIZE(MasterEntities)
1836
str = TRIM(str)//' '//I2S(MasterEntities(i))
1837
END DO
1838
CALL Info(Caller,'Performing reduction for BCs: '//TRIM(str),Level=10)
1839
1840
Omega = ListGetAngularFrequency(Found=Found)
1841
IF(.NOT. Found) CALL Fatal(Caller,'We need angular frequency!')
1842
1843
t1 = Mesh % NumberOfBulkElements + 1
1844
t2 = Mesh % NumberOfBulkElements + Mesh % NumberOfBoundaryElements
1845
1846
EdgeBasis = .FALSE.
1847
UseGaussLaw = .FALSE.
1848
1849
IF(avar % dofs <= 2) THEN
1850
EdgeBasis = .TRUE.
1851
CALL EdgeElementStyle(avar % Solver % Values, PiolaVersion, BasisDegree = EdgeBasisDegree )
1852
UseGaussLaw = ListGetLogical(avar % solver % values, 'Use Gauss Law', Found)
1853
END IF
1854
1855
OutFlux = 0.0_dp
1856
area = 0.0_dp
1857
1858
IF( ALL([1,2,3,6] /= Avar % Dofs) ) THEN
1859
CALL Fatal(Caller,'Invalid number of components for vector potential!')
1860
END IF
1861
1862
InitHandles = .TRUE.
1863
DO t=t1, t2
1864
Element => Mesh % Elements(t)
1865
bc_id = Element % BoundaryInfo % Constraint
1866
IF( ALL( MasterEntities /= bc_id ) ) CYCLE
1867
BC => Model % BCs(bc_id) % Values
1868
1869
IF(EdgeBasis .AND. Element % Type % ElementCode > 300 ) THEN
1870
k = FindBoundaryFaceIndex(Mesh,Element)
1871
Element => Mesh % Faces(k)
1872
END IF
1873
IF(UseGaussLaw) THEN
1874
CALL LocalIntegBC_AV(BC, Element, InitHandles)
1875
ELSE
1876
CALL LocalIntegBC_E(BC, Element, InitHandles)
1877
END IF
1878
END DO
1879
1880
Area = ParallelReduction(Area)
1881
trans = ParallelReduction(trans)
1882
Zimp = 1.0_dp / trans
1883
1884
PortImp = Zimp
1885
1886
IF( UseGaussLaw ) THEN
1887
vol = ParallelReduction(Vol)
1888
curr = ParallelReduction(curr)
1889
port_curr = ParallelReduction(port_curr)
1890
1891
! Still testing these!
1892
curr = -curr
1893
port_curr = -port_curr
1894
1895
vol = vol / area
1896
1897
PRINT *,'LumpedCurr av:',vol,curr,port_curr,Zimp *curr, CONJG(Zimp) * port_curr
1898
1899
OutFlux = (vol + Zimp * curr) / (2*SQRT(REAL(Zimp)))
1900
InFlux = (vol - CONJG(Zimp) * port_curr ) / (2*SQRT(REAL(Zimp)))
1901
1902
! For now use just the average voltages
1903
OutFlux = vol !curr
1904
InFlux = 1.0 !port_curr
1905
PortImp = Zimp
1906
ELSE
1907
int_el = ParallelReduction(int_el)
1908
int_norm = ParallelReduction(int_norm)
1909
OutFlux = int_el
1910
InFlux = int_norm
1911
1912
IF(ABS(int_norm) < 1.0e-20 ) THEN
1913
PRINT *,'int_norm:',int_norm,area,EPSILON(area)
1914
CALL Warn(Caller,'Source seems to be close to zero!')
1915
END IF
1916
1917
PRINT *,'LumpedCurr e:',int_el,int_norm,area,trans,Zimp,ABS(int_norm)
1918
1919
END IF
1920
1921
1922
CALL Info(Caller,'Reduction operator finished',Level=12)
1923
1924
CONTAINS
1925
1926
!-----------------------------------------------------------------------------
1927
SUBROUTINE LocalIntegBC_E( BC, Element, InitHandles )
1928
!------------------------------------------------------------------------------
1929
TYPE(ValueList_t), POINTER :: BC
1930
TYPE(Element_t), POINTER :: Element
1931
LOGICAL :: InitHandles
1932
!------------------------------------------------------------------------------
1933
COMPLEX(KIND=dp) :: B, L(3), muinv, MagLoad(3), TemGrad(3), &
1934
e_ip(3), e_ip_norm, e_ip_tan(3), f_ip_tan(3), imu, phi, eps0, mu0inv, epsr, mur, ElSurfCurr(3)
1935
REAL(KIND=dp), ALLOCATABLE :: Basis(:),dBasisdx(:,:),WBasis(:,:),RotWBasis(:,:), e_local(:,:)
1936
REAL(KIND=dp) :: weight, DetJ, Normal(3), cond, u, v, w, x, y, z, rob0
1937
TYPE(Nodes_t), SAVE :: ElementNodes, ParentNodes
1938
INTEGER, POINTER :: NodeIndexes(:), pIndexes(:), ParentIndexes(:)
1939
INTEGER, POINTER, SAVE :: EdgeIndexes(:)
1940
LOGICAL :: Stat, Found
1941
TYPE(GaussIntegrationPoints_t) :: IP
1942
INTEGER :: t, i, j, m, np, p, q, ndofs, n, nd
1943
LOGICAL :: AllocationsDone = .FALSE.
1944
TYPE(Element_t), POINTER :: Parent, Parent2
1945
TYPE(ValueHandle_t), SAVE :: MagLoad_h, ElRobin_h, MuCoeff_h, Absorb_h, TemRe_h, TemIm_h, ElSurfCurr_h
1946
TYPE(ValueHandle_t), SAVE :: CondCoeff_h, CurrDens_h, EpsCoeff_h
1947
INTEGER :: nactive
1948
1949
TYPE(ValueHandle_t), SAVE :: PortTypeIndex_h, PortZ_h, PortLength_h, PortScale_h, PortDirection_h, PortCenter_h
1950
INTEGER :: PortTypeIndex, PortDirection
1951
COMPLEX(KIND=dp) :: PortZ
1952
REAL(KIND=dp) :: PortLength, PortScale, PortCenter(3)
1953
LOGICAL :: GotPort
1954
1955
SAVE AllocationsDone, WBasis, RotWBasis, Basis, dBasisdx, e_local, mu0inv, eps0
1956
1957
ndofs = avar % dofs
1958
IF(.NOT. AllocationsDone ) THEN
1959
m = Mesh % MaxElementDOFs
1960
ALLOCATE( ElementNodes % x(m), ElementNodes % y(m), ElementNodes % z(m), EdgeIndexes(m), &
1961
ParentNodes % x(m), ParentNodes % y(m), ParentNodes % z(m), &
1962
WBasis(m,3), RotWBasis(m,3), Basis(m), dBasisdx(m,3), e_local(ndofs,m) )
1963
AllocationsDone = .TRUE.
1964
END IF
1965
1966
IF( InitHandles ) THEN
1967
CALL ListInitElementKeyword( ElRobin_h,'Boundary Condition','Electric Robin Coefficient',InitIm=.TRUE.)
1968
CALL ListInitElementKeyword( MagLoad_h,'Boundary Condition','Magnetic Boundary Load', InitIm=.TRUE.,InitVec3D=.TRUE.)
1969
CALL ListInitElementKeyword( ElSurfCurr_h, 'Boundary Condition', 'Electric Surface Current', &
1970
InitIm = .TRUE., InitVec3D=.TRUE.)
1971
CALL ListInitElementKeyword( Absorb_h,'Boundary Condition','Absorbing BC')
1972
CALL ListInitElementKeyword( TemRe_h,'Boundary Condition','TEM Potential')
1973
CALL ListInitElementKeyword( TemIm_h,'Boundary Condition','TEM Potential Im')
1974
1975
CALL ListInitElementKeyword( MuCoeff_h,'Material','Relative Reluctivity',InitIm=.TRUE.)
1976
CALL ListInitElementKeyword( EpsCoeff_h,'Material','Relative Permittivity',InitIm=.TRUE.)
1977
CALL ListInitElementKeyword( CondCoeff_h,'Material','Electric Conductivity')
1978
1979
! Lumped ports
1980
CALL ListInitElementKeyword( PortTypeIndex_h,'Boundary Condition','Port Type Index')
1981
CALL ListInitElementKeyword( PortZ_h,'Boundary Condition','Port Impedance',InitIm=.TRUE.)
1982
CALL ListInitElementKeyword( PortLength_h,'Boundary Condition','Port Length')
1983
CALL ListInitElementKeyword( PortScale_h,'Boundary Condition','Port Scale')
1984
CALL ListInitElementKeyword( PortDirection_h,'Boundary Condition','Port Direction',DefIValue=3)
1985
CALL ListInitElementKeyword( PortCenter_h,'Boundary Condition','Port Center',InitVec3D=.TRUE.)
1986
1987
Found = .FALSE.
1988
IF( ASSOCIATED( Model % Constants ) ) THEN
1989
mu0inv = ListGetConstReal( Model % Constants,'Permeability of Vacuum', Found )
1990
IF(mu0inv/=0) mu0inv = 1/mu0inv;
1991
END IF
1992
IF(.NOT. Found ) mu0inv = 1.0_dp / ( PI * 4.0d-7 )
1993
Found = .FALSE.
1994
IF( ASSOCIATED( Model % Constants ) ) THEN
1995
eps0 = ListGetConstReal ( Model % Constants,'Permittivity of Vacuum', Found )
1996
END IF
1997
IF(.NOT. Found ) eps0 = 8.854187817d-12
1998
InitHandles = .FALSE.
1999
END IF
2000
2001
imu = CMPLX(0.0_dp, 1.0_dp, KIND=dp)
2002
rob0 = Omega * SQRT( eps0 / mu0inv )
2003
2004
n = Element % TYPE % NumberOfNodes
2005
NodeIndexes => Element % NodeIndexes
2006
2007
ElementNodes % x(1:n) = Mesh % Nodes % x(NodeIndexes(1:n))
2008
ElementNodes % y(1:n) = Mesh % Nodes % y(NodeIndexes(1:n))
2009
ElementNodes % z(1:n) = Mesh % Nodes % z(NodeIndexes(1:n))
2010
2011
! We always need the parent for material properties.
2012
Parent => Element % BoundaryInfo % Left
2013
IF(.NOT. ASSOCIATED(Parent)) THEN
2014
Parent => Element % BoundaryInfo % Right
2015
END IF
2016
IF(.NOT. ASSOCIATED( Parent ) ) THEN
2017
CALL Fatal(Caller,'Model lumping requires parent element!')
2018
END IF
2019
2020
! If the source of the incident field is not in the plane of the port we have to use the parent!
2021
#define doparent 0
2022
IF( EdgeBasis ) THEN
2023
#if doparent
2024
np = Parent % TYPE % NumberOfNodes
2025
ParentIndexes => Parent % NodeIndexes
2026
ParentNodes % x(1:np) = Mesh % Nodes % x(ParentIndexes(1:np))
2027
ParentNodes % y(1:np) = Mesh % Nodes % y(ParentIndexes(1:np))
2028
ParentNodes % z(1:np) = Mesh % Nodes % z(ParentIndexes(1:np))
2029
2030
nd = mGetElementDofs( EdgeIndexes, Uelement = Parent, USolver = avar % Solver )
2031
np = COUNT(EdgeIndexes(1:nd) <= Mesh % NumberOfNodes)
2032
pIndexes => EdgeIndexes
2033
#else
2034
nd = mGetElementDofs( EdgeIndexes, Uelement = Element, USolver = avar % Solver )
2035
np = COUNT(EdgeIndexes(1:nd) <= Mesh % NumberOfNodes)
2036
pIndexes => EdgeIndexes
2037
2038
#endif
2039
IP = GaussPoints(Element, EdgeBasis=.TRUE., PReferenceElement=PiolaVersion, &
2040
EdgeBasisDegree=EdgeBasisDegree)
2041
ELSE
2042
nd = n
2043
IP = GaussPoints( Element )
2044
IF( avar % TYPE == Variable_on_nodes_on_elements ) THEN
2045
pIndexes => Element % DGIndexes
2046
ELSE
2047
pIndexes => NodeIndexes
2048
END IF
2049
END IF
2050
2051
! Edge or nodal values of vector potential
2052
DO i=1,avar % dofs
2053
e_local(i,1:nd) = avar % values(avar % dofs*(avar % Perm(pIndexes(1:nd))-1)+i)
2054
END DO
2055
2056
PortTypeIndex = ListGetElementInteger(PortTypeIndex_h, Element, GotPort)
2057
IF(GotPort) THEN
2058
PortTypeIndex = ListGetElementInteger( PortTypeIndex_h, Element )
2059
PortZ = ListGetElementComplex( PortZ_h, Element = Element )
2060
PortScale = ListGetElementReal( PortScale_h, Element = Element )
2061
PortLength = ListGetElementReal( PortLength_h, Element = Element )
2062
IF( PortTypeIndex == 1 ) THEN
2063
PortDirection = ListGetElementInteger( PortDirection_h, Element )
2064
ELSE
2065
PortCenter = ListGetElementReal( PortCenter_h, Element = Element )
2066
END IF
2067
END IF
2068
2069
#if doparent
2070
! Normal is only needed if we integrate over parent element.
2071
Normal = NormalVector(Element, ElementNodes, Check=.TRUE.)
2072
#endif
2073
2074
! Numerical integration:
2075
!-----------------------
2076
DO t=1,IP % n
2077
2078
stat = ElementInfo( Element, ElementNodes, IP % U(t), IP % V(t), &
2079
IP % W(t), detJ, Basis, dBasisdx )
2080
weight = IP % s(t) * detJ
2081
2082
! Get material properties from parent element.
2083
!----------------------------------------------
2084
mur = ListGetElementComplex( MuCoeff_h, Basis, Parent, Found, GaussPoint = t )
2085
IF( .NOT. Found ) mur = 1.0_dp
2086
muinv = mur * mu0inv
2087
2088
Cond = ListGetElementReal( CondCoeff_h, Basis, Parent, Found, GaussPoint = t )
2089
2090
L = (0_dp, 0_dp)
2091
2092
IF( ListGetElementLogical( Absorb_h, Element, Found ) ) THEN
2093
epsr = ListGetElementComplex( EpsCoeff_h, Basis, Parent, Found, GaussPoint = t )
2094
IF( .NOT. Found ) epsr = 1.0_dp
2095
B = imu * rob0 * SQRT( epsr / mur )
2096
2097
ELSE IF(GotPort) THEN
2098
IF( PortTypeIndex == 1 ) THEN
2099
B = CMPLX(0_dp, 1_dp, KIND=dp) * ( omega / mu0inv ) / (PortScale * PortZ )
2100
L(ABS(PortDirection)) = SIGN(1,PortDirection) / ( PortLength * SQRT( PortScale ) )
2101
END IF
2102
ELSE
2103
B = ListGetElementComplex( ElRobin_h, Basis, Element, Found, GaussPoint = t )
2104
2105
MagLoad = ListGetElementComplex3D( MagLoad_h, Basis, Element, Found, GaussPoint = t )
2106
ElSurfCurr = ListGetElementComplex3D( ElSurfCurr_h, Basis, Element, Found, GaussPoint = t)
2107
2108
TemGrad = CMPLX( ListGetElementRealGrad( TemRe_h,dBasisdx,Element,Found), &
2109
ListGetElementRealGrad( TemIm_h,dBasisdx,Element,Found), KIND=dp )
2110
2111
IF (ABS(B) > AEPS) THEN
2112
L = ( MagLoad + TemGrad ) / ( 2*B)
2113
END IF
2114
2115
! We sum the components here, since continutation is a little cleaner then...
2116
ElSurfCurr = L + (0_dp, 1_dp)*omega*ElSurfCurr
2117
END IF
2118
2119
2120
IF( EdgeBasis ) THEN
2121
#if doparent
2122
! In order to get the normal component of the electric field we must operate on the
2123
! parent element. The surface element only has tangential components.
2124
! Note: the finding of parents does not work for piola/quadratic gauss points.
2125
CALL FindParentUVW( Element, n, Parent, Parent % TYPE % NumberOfNodes, U, V, W, Basis )
2126
stat = ElementInfo( Parent, ParentNodes, u, v, w, detJ, Basis, dBasisdx, &
2127
EdgeBasis = Wbasis, RotBasis = RotWBasis, USolver = avar % Solver )
2128
#else
2129
stat = ElementInfo( Element, ElementNodes, IP % U(t), IP % V(t), &
2130
IP % W(t), detJ, Basis, dBasisdx, &
2131
EdgeBasis = Wbasis, RotBasis = RotWBasis, USolver = avar % Solver )
2132
#endif
2133
e_ip(1:3) = CMPLX(MATMUL(e_local(1,np+1:nd),WBasis(1:nd-np,1:3)), MATMUL(e_local(2,np+1:nd),WBasis(1:nd-np,1:3)), KIND=dp)
2134
ELSE
2135
DO i=1,3
2136
e_ip(i) = CMPLX( SUM( Basis(1:n) * e_local(i,1:n) ), SUM( Basis(1:n) * e_local(i+3,1:n) ), KIND=dp )
2137
END DO
2138
END IF
2139
2140
! Integral over electric field: This gives the phase
2141
#if doparent
2142
! e_ip is actually by construction on the plane if we do not use parent elements
2143
e_ip_norm = SUM(e_ip*Normal)
2144
e_ip_tan = e_ip - e_ip_norm * Normal
2145
int_el = int_el + weight * SUM(e_ip_tan * CONJG(L) )
2146
#else
2147
int_el = int_el + weight * SUM(e_ip * CONJG(L) )
2148
#endif
2149
2150
! Norm of electric field used for normalization
2151
int_norm = int_norm + weight * ABS( SUM( L * CONJG(L) ) )
2152
2153
trans = trans + B * weight / Omega
2154
area = area + weight
2155
END DO
2156
2157
!------------------------------------------------------------------------------
2158
END SUBROUTINE LocalIntegBC_E
2159
!------------------------------------------------------------------------------
2160
2161
FUNCTION RealComplexCrossProduct(v1,v2) RESULT(v3)
2162
!------------------------------------------------------------------------------
2163
COMPLEX(KIND=dp) :: v2(3), v3(3)
2164
REAL(KIND=dp) :: v1(3)
2165
v3(1) = v1(2)*v2(3) - v1(3)*v2(2)
2166
v3(2) = -v1(1)*v2(3) + v1(3)*v2(1)
2167
v3(3) = v1(1)*v2(2) - v1(2)*v2(1)
2168
!------------------------------------------------------------------------------
2169
END FUNCTION RealComplexCrossProduct
2170
2171
!-----------------------------------------------------------------------------
2172
SUBROUTINE LocalIntegBC_AV( BC, Element, InitHandles )
2173
!------------------------------------------------------------------------------
2174
TYPE(ValueList_t), POINTER :: BC
2175
TYPE(Element_t), POINTER :: Element
2176
LOGICAL :: InitHandles
2177
!------------------------------------------------------------------------------
2178
COMPLEX(KIND=dp) :: tc_ip, cd_ip, v_ip, ep_ip, mu0inv, muinv, mur, &
2179
cond_ip, imu
2180
REAL(KIND=dp), ALLOCATABLE :: Basis(:),dBasisdx(:,:),v_local(:,:)
2181
REAL(KIND=dp) :: weight, DetJ
2182
TYPE(Nodes_t), SAVE :: ElementNodes
2183
INTEGER, POINTER :: NodeIndexes(:), pIndexes(:), ParentIndexes(:)
2184
LOGICAL :: Found
2185
TYPE(GaussIntegrationPoints_t) :: IP
2186
INTEGER :: t, i, j, m, np, p, q, ndofs, n, nd
2187
LOGICAL :: AllocationsDone = .FALSE.
2188
TYPE(Element_t), POINTER :: Parent, MatElement
2189
TYPE(ValueHandle_t), SAVE :: MuCoeff_h, CondCoeff_h, ExtPot_h
2190
TYPE(ValueHandle_t), SAVE :: TransferCoeff_h, ElCurrent_h, BCMat_h
2191
2192
SAVE AllocationsDone, Basis, dBasisdx, v_local, mu0inv
2193
2194
ndofs = avar % dofs
2195
IF(.NOT. AllocationsDone ) THEN
2196
m = Mesh % MaxElementDOFs
2197
ALLOCATE( ElementNodes % x(m), ElementNodes % y(m), ElementNodes % z(m), &
2198
Basis(m), dBasisdx(m,3), v_local(ndofs,m) )
2199
AllocationsDone = .TRUE.
2200
END IF
2201
2202
! BC given with these:
2203
! Electric Transfer Coefficient
2204
! Electric Current Density / Incident Voltage
2205
IF( InitHandles ) THEN
2206
CALL ListInitElementKeyword( MuCoeff_h,'Material','Relative Reluctivity',InitIm=.TRUE.)
2207
CALL ListInitElementKeyword( CondCoeff_h,'Material','Electric Conductivity')
2208
2209
CALL ListInitElementKeyword( TransferCoeff_h,'Boundary Condition','Electric Transfer Coefficient',InitIm=.TRUE.)
2210
CALL ListInitElementKeyword( ElCurrent_h,'Boundary Condition','Electric Current Density',InitIm=.TRUE.)
2211
CALL ListInitElementKeyword( ExtPot_h,'Boundary Condition','Incident Voltage',InitIm=.TRUE.)
2212
2213
CALL ListInitElementKeyword( BCMat_h,'Boundary Condition','Material')
2214
2215
Found = .FALSE.
2216
IF( ASSOCIATED( Model % Constants ) ) THEN
2217
mu0inv = 1.0_dp / ListGetConstReal( Model % Constants,'Permeability of Vacuum', Found )
2218
END IF
2219
IF(.NOT. Found ) mu0inv = 1.0_dp / ( PI * 4.0d-7 )
2220
Found = .FALSE.
2221
InitHandles = .FALSE.
2222
END IF
2223
2224
imu = CMPLX(0.0_dp, 1.0_dp, KIND=dp)
2225
2226
n = Element % TYPE % NumberOfNodes
2227
NodeIndexes => Element % NodeIndexes
2228
2229
ElementNodes % x(1:n) = Mesh % Nodes % x(NodeIndexes(1:n))
2230
ElementNodes % y(1:n) = Mesh % Nodes % y(NodeIndexes(1:n))
2231
ElementNodes % z(1:n) = Mesh % Nodes % z(NodeIndexes(1:n))
2232
2233
Parent => Element % BoundaryInfo % Left
2234
IF(.NOT. ASSOCIATED(Parent)) Parent => Element % BoundaryInfo % Right
2235
IF(.NOT. ASSOCIATED( Parent ) ) THEN
2236
CALL Fatal(Caller,'Model lumping requires parent element!')
2237
END IF
2238
2239
! If we have material also define in the BC section then use it
2240
i = ListGetElementInteger( BCMat_h, Element, Found )
2241
IF( i > 0 ) THEN
2242
MatElement => Element
2243
ELSE
2244
MatElement => Parent
2245
END IF
2246
2247
nd = n
2248
IP = GaussPoints( Element )
2249
pIndexes => NodeIndexes
2250
2251
! Nodal values of scalar potential
2252
DO i=1,avar % dofs
2253
v_local(i,1:nd) = avar % values(avar % dofs*(avar % Perm(pIndexes(1:nd))-1)+i)
2254
END DO
2255
2256
! Numerical integration:
2257
!-----------------------
2258
DO t=1,IP % n
2259
2260
stat = ElementInfo( Element, ElementNodes, IP % U(t), IP % V(t), &
2261
IP % W(t), detJ, Basis, dBasisdx )
2262
weight = IP % s(t) * detJ
2263
2264
! Get material properties from parent element.
2265
!----------------------------------------------
2266
mur = ListGetElementComplex( MuCoeff_h, Basis, MatElement, Found, GaussPoint = t )
2267
IF( .NOT. Found ) mur = 1.0_dp
2268
muinv = mur * mu0inv
2269
2270
cond_ip = ListGetElementReal( CondCoeff_h, Basis, MatElement, Found, GaussPoint = t )
2271
cd_ip = ListGetElementComplex( ElCurrent_h, Basis, Element, Found, GaussPoint = t )
2272
2273
tc_ip = ListGetElementComplex( TransferCoeff_h, Basis, Element, Found, GaussPoint = t )
2274
IF(Found) THEN
2275
ep_ip = ListGetElementComplex( ExtPot_h, Basis, Element, Found, GaussPoint = t )
2276
IF(Found) cd_ip = cd_ip + 2 * tc_ip * ep_ip
2277
END IF
2278
v_ip = CMPLX( SUM( Basis(1:n) * v_local(1,1:n) ), SUM( Basis(1:n) * v_local(2,1:n) ), KIND=dp )
2279
2280
area = area + weight
2281
2282
curr = curr - tc_ip * v_ip * weight
2283
port_curr = port_curr + cd_ip * weight
2284
trans = trans + tc_ip * weight
2285
vol = vol + v_ip * weight
2286
END DO
2287
2288
!------------------------------------------------------------------------------
2289
END SUBROUTINE LocalIntegBC_AV
2290
!------------------------------------------------------------------------------
2291
2292
END FUNCTION BoundaryWaveFlux
2293
2294
2295
2296
!------------------------------------------------------------------------------
2297
!> This routine computes geometric quantities for lumped port making the setting
2298
!> of BCs a little easier.
2299
!------------------------------------------------------------------------------
2300
SUBROUTINE DefinePortParameters(Model, Mesh)
2301
!------------------------------------------------------------------------------
2302
IMPLICIT NONE
2303
TYPE(Model_t) :: Model
2304
TYPE(Mesh_t), POINTER :: Mesh
2305
!------------------------------------------------------------------------------
2306
TYPE(Nodes_t) :: ElementNodes
2307
TYPE(Element_t), POINTER :: Element
2308
TYPE(GaussIntegrationPoints_t) :: IP
2309
INTEGER ::i,j,k,n,t,t0,bc_id,ierr,PortDir,PortIndex,PortTypeInd
2310
INTEGER, POINTER :: Indexes(:)
2311
REAL(KIND=dp) :: s,detJ,Width,Area,Length,RadInner, RadOuter,CenterArray(3,1),Scale
2312
REAL(KIND=dp), ALLOCATABLE :: Basis(:), LumpVec(:)
2313
TYPE(ValueList_t), POINTER :: BC, BC0
2314
CHARACTER(:), ALLOCATABLE :: PortType
2315
LOGICAL :: Found, stat
2316
2317
CHARACTER(*), PARAMETER :: Caller = 'DefinePortParameters'
2318
2319
IF(.NOT. ASSOCIATED( Mesh ) ) THEN
2320
CALL Fatal(Caller,'Mesh not associated!')
2321
END IF
2322
2323
CALL Info(Caller,'Defining geometric port parameters',Level=8)
2324
2325
n = Mesh % MaxElementNodes
2326
t0 = Mesh % NumberOfBulkElements
2327
2328
! Check the number of ports and add "port type index" keyword.
2329
DO bc_id = 1,Model % NumberOfBCs
2330
BC => Model % BCs(bc_id) % Values
2331
2332
! This has already been defined!
2333
IF( ListCheckPresent( BC,'Port Type Index')) CYCLE
2334
2335
IF(.NOT. ListCheckPresent(BC,'port impedance')) CYCLE
2336
2337
PortType = ListGetString( BC,'port type',Found)
2338
2339
SELECT CASE(PortType)
2340
CASE('rectangular')
2341
PortTypeInd = 1
2342
2343
CASE('coaxial')
2344
PortTypeInd = 2
2345
2346
CASE DEFAULT
2347
CALL Info(Caller,'Port Type "Port Type" defaulted to "rectangular"',Level=4)
2348
PortTypeInd = 1
2349
END SELECT
2350
2351
CALL Info(Caller,'Defining parameters for port on BC: '//I2S(bc_id),Level=8)
2352
CALL ListAddInteger( BC,'Port Type Index',PortTypeInd)
2353
2354
IF(.NOT. ALLOCATED(Basis) ) THEN
2355
ALLOCATE(Basis(n), ElementNodes % x(n), ElementNodes % y(n), ElementNodes % z(n), LumpVec(7) )
2356
END IF
2357
2358
LumpVec(1:3) = HUGE(s)
2359
LumpVec(4:6) = -HUGE(s)
2360
LumpVec(7) = 0.0_dp
2361
2362
2363
DO t=1, Mesh % NumberOfBoundaryElements
2364
Element => Mesh % Elements( t0 + t )
2365
IF( Element % BoundaryInfo % Constraint /= Model % BCs(bc_id) % Tag ) CYCLE
2366
2367
Indexes => Element % NodeIndexes
2368
n = Element % TYPE % NumberOfNodes
2369
2370
ElementNodes % x = 0.0_dp
2371
ElementNodes % y = 0.0_dp
2372
ElementNodes % z = 0.0_dp
2373
2374
ElementNodes % x(1:n) = Mesh % Nodes % x(Indexes)
2375
ElementNodes % y(1:n) = Mesh % Nodes % y(Indexes)
2376
ElementNodes % z(1:n) = Mesh % Nodes % z(Indexes)
2377
2378
! Get min/max range for each coordinate.
2379
LumpVec(1) = MIN(LumpVec(1),MINVAL(ElementNodes % x(1:n)))
2380
LumpVec(2) = MIN(LumpVec(2),MINVAL(ElementNodes % y(1:n)))
2381
LumpVec(3) = MIN(LumpVec(3),MINVAL(ElementNodes % z(1:n)))
2382
LumpVec(4) = MAX(LumpVec(4),MAXVAL(ElementNodes % x(1:n)))
2383
LumpVec(5) = MAX(LumpVec(5),MAXVAL(ElementNodes % y(1:n)))
2384
LumpVec(6) = MAX(LumpVec(6),MAXVAL(ElementNodes % z(1:n)))
2385
2386
! Integrate over the area.
2387
IP = GaussPoints( Element, PReferenceElement = .FALSE.)
2388
DO j=1,IP % n
2389
stat = ElementInfo( Element, ElementNodes, IP % U(j), IP % V(j), IP % W(j), detJ, Basis )
2390
S = DetJ * IP % s(j)
2391
LumpVec(7) = LumpVec(7) + s
2392
END DO
2393
END DO
2394
2395
! Do parallel communication, if needed.
2396
IF( ParEnv % PEs > 1 ) THEN
2397
CALL MPI_ALLREDUCE( MPI_IN_PLACE, LumpVec(1:3), 3, &
2398
MPI_DOUBLE_PRECISION, MPI_MIN, ELMER_COMM_WORLD, ierr )
2399
CALL MPI_ALLREDUCE( MPI_IN_PLACE, LumpVec(4:6), 3, &
2400
MPI_DOUBLE_PRECISION, MPI_MAX, ELMER_COMM_WORLD, ierr )
2401
CALL MPI_ALLREDUCE( MPI_IN_PLACE, LumpVec(7:7), 1, &
2402
MPI_DOUBLE_PRECISION, MPI_SUM, ELMER_COMM_WORLD, ierr )
2403
END IF
2404
2405
! Area is used in all port models.
2406
Area = LumpVec(7)
2407
2408
SELECT CASE(PortTypeInd)
2409
CASE(1)
2410
PortDir = ABS( ListGetInteger( BC,'Port Direction',Found) )
2411
IF(.NOT. Found) PortDir = 3
2412
2413
Length = LumpVec(3+PortDir) - LumpVec(PortDir)
2414
Width = Area / Length
2415
Scale = Width / Length
2416
2417
!PRINT *,'area:',area, length, width, scale
2418
2419
CALL ListAddConstReal( BC,'Port Length',Length)
2420
CALL ListAddConstReal( BC,'Port Scale',Scale)
2421
IF(InfoActive(8)) THEN
2422
PRINT *,'Setting rectangular port parameters:',Length,Scale
2423
END IF
2424
2425
CASE(2)
2426
RadOuter = 0.0_dp
2427
DO i=1,3
2428
RadOuter = MAX(RadOuter,(LumpVec(3+i)-LumpVec(i))/2)
2429
CenterArray(i,1) = (LumpVec(3+i)+LumpVec(i))/2
2430
END DO
2431
RadInner = SQRT(RadOuter**2-Area/PI)
2432
Length = (RadInner+RadOuter)/2
2433
Scale = 2*PI/LOG(RadOuter/RadInner)
2434
2435
! PRINT *,'area:',area, radinner, radouter
2436
2437
CALL ListAddConstReal( BC,'Port Length',Length)
2438
CALL ListAddConstReal( BC,'Port Scale',Scale)
2439
CALL ListAddConstRealArray( BC,'Port Center',3,1,CenterArray)
2440
IF(InfoActive(8)) THEN
2441
PRINT *,'Setting coaxial port parameters:',Length,Scale,' and center ',CenterArray
2442
END IF
2443
END SELECT
2444
2445
END DO
2446
2447
IF(ALLOCATED(Basis)) THEN
2448
DEALLOCATE(Basis, ElementNodes % x, ElementNodes % y, ElementNodes % z, LumpVec )
2449
END IF
2450
2451
END SUBROUTINE DefinePortParameters
2452
!------------------------------------------------------------------------------
2453
2454
2455
2456
2457
END MODULE LumpingUtils
2458
!------------------------------------------------------------------------------
2459
2460
2461
!> \}
2462
2463
2464
2465