-
Notifications
You must be signed in to change notification settings - Fork 83
/
lulesh.cc
2792 lines (2332 loc) · 89.1 KB
/
lulesh.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
This is a Version 2.0 MPI + OpenMP implementation of LULESH
Copyright (c) 2010-2013.
Lawrence Livermore National Security, LLC.
Produced at the Lawrence Livermore National Laboratory.
LLNL-CODE-461231
All rights reserved.
This file is part of LULESH, Version 2.0.
Please also read this link -- http://www.opensource.org/licenses/index.php
//////////////
DIFFERENCES BETWEEN THIS VERSION (2.x) AND EARLIER VERSIONS:
* Addition of regions to make work more representative of multi-material codes
* Default size of each domain is 30^3 (27000 elem) instead of 45^3. This is
more representative of our actual working set sizes
* Single source distribution supports pure serial, pure OpenMP, MPI-only,
and MPI+OpenMP
* Addition of ability to visualize the mesh using VisIt
https://wci.llnl.gov/codes/visit/download.html
* Various command line options (see ./lulesh2.0 -h)
-q : quiet mode - suppress stdout
-i <iterations> : number of cycles to run
-s <size> : length of cube mesh along side
-r <numregions> : Number of distinct regions (def: 11)
-b <balance> : Load balance between regions of a domain (def: 1)
-c <cost> : Extra cost of more expensive regions (def: 1)
-f <filepieces> : Number of file parts for viz output (def: np/9)
-p : Print out progress
-v : Output viz file (requires compiling with -DVIZ_MESH
-h : This message
printf("Usage: %s [opts]\n", execname);
printf(" where [opts] is one or more of:\n");
printf(" -q : quiet mode - suppress all stdout\n");
printf(" -i <iterations> : number of cycles to run\n");
printf(" -s <size> : length of cube mesh along side\n");
printf(" -r <numregions> : Number of distinct regions (def: 11)\n");
printf(" -b <balance> : Load balance between regions of a domain (def: 1)\n");
printf(" -c <cost> : Extra cost of more expensive regions (def: 1)\n");
printf(" -f <numfiles> : Number of files to split viz dump into (def: (np+10)/9)\n");
printf(" -p : Print out progress\n");
printf(" -v : Output viz file (requires compiling with -DVIZ_MESH\n");
printf(" -h : This message\n");
printf("\n\n");
*Notable changes in LULESH 2.0
* Split functionality into different files
lulesh.cc - where most (all?) of the timed functionality lies
lulesh-comm.cc - MPI functionality
lulesh-init.cc - Setup code
lulesh-viz.cc - Support for visualization option
lulesh-util.cc - Non-timed functions
*
* The concept of "regions" was added, although every region is the same ideal
* gas material, and the same sedov blast wave problem is still the only
* problem its hardcoded to solve.
* Regions allow two things important to making this proxy app more representative:
* Four of the LULESH routines are now performed on a region-by-region basis,
* making the memory access patterns non-unit stride
* Artificial load imbalances can be easily introduced that could impact
* parallelization strategies.
* The load balance flag changes region assignment. Region number is raised to
* the power entered for assignment probability. Most likely regions changes
* with MPI process id.
* The cost flag raises the cost of ~45% of the regions to evaluate EOS by the
* entered multiple. The cost of 5% is 10x the entered multiple.
* MPI and OpenMP were added, and coalesced into a single version of the source
* that can support serial builds, MPI-only, OpenMP-only, and MPI+OpenMP
* Added support to write plot files using "poor mans parallel I/O" when linked
* with the silo library, which in turn can be read by VisIt.
* Enabled variable timestep calculation by default (courant condition), which
* results in an additional reduction.
* Default domain (mesh) size reduced from 45^3 to 30^3
* Command line options to allow numerous test cases without needing to recompile
* Performance optimizations and code cleanup beyond LULESH 1.0
* Added a "Figure of Merit" calculation (elements solved per microsecond) and
* output in support of using LULESH 2.0 for the 2017 CORAL procurement
*
* Possible Differences in Final Release (other changes possible)
*
* High Level mesh structure to allow data structure transformations
* Different default parameters
* Minor code performance changes and cleanup
TODO in future versions
* Add reader for (truly) unstructured meshes, probably serial only
* CMake based build system
//////////////
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions
are met:
* Redistributions of source code must retain the above copyright
notice, this list of conditions and the disclaimer below.
* Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the disclaimer (as noted below)
in the documentation and/or other materials provided with the
distribution.
* Neither the name of the LLNS/LLNL nor the names of its contributors
may be used to endorse or promote products derived from this software
without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL LAWRENCE LIVERMORE NATIONAL SECURITY, LLC,
THE U.S. DEPARTMENT OF ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT,
INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY
OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE,
EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
Additional BSD Notice
1. This notice is required to be provided under our contract with the U.S.
Department of Energy (DOE). This work was produced at Lawrence Livermore
National Laboratory under Contract No. DE-AC52-07NA27344 with the DOE.
2. Neither the United States Government nor Lawrence Livermore National
Security, LLC nor any of their employees, makes any warranty, express
or implied, or assumes any liability or responsibility for the accuracy,
completeness, or usefulness of any information, apparatus, product, or
process disclosed, or represents that its use would not infringe
privately-owned rights.
3. Also, reference herein to any specific commercial products, process, or
services by trade name, trademark, manufacturer or otherwise does not
necessarily constitute or imply its endorsement, recommendation, or
favoring by the United States Government or Lawrence Livermore National
Security, LLC. The views and opinions of authors expressed herein do not
necessarily state or reflect those of the United States Government or
Lawrence Livermore National Security, LLC, and shall not be used for
advertising or product endorsement purposes.
*/
#include <climits>
#include <vector>
#include <math.h>
#include <stdio.h>
#include <string.h>
#include <ctype.h>
#include <time.h>
#include <sys/time.h>
#include <iostream>
#include <unistd.h>
#if _OPENMP
# include <omp.h>
#endif
#include "lulesh.h"
/* Work Routines */
static inline
void TimeIncrement(Domain& domain)
{
Real_t targetdt = domain.stoptime() - domain.time() ;
if ((domain.dtfixed() <= Real_t(0.0)) && (domain.cycle() != Int_t(0))) {
Real_t ratio ;
Real_t olddt = domain.deltatime() ;
/* This will require a reduction in parallel */
Real_t gnewdt = Real_t(1.0e+20) ;
Real_t newdt ;
if (domain.dtcourant() < gnewdt) {
gnewdt = domain.dtcourant() / Real_t(2.0) ;
}
if (domain.dthydro() < gnewdt) {
gnewdt = domain.dthydro() * Real_t(2.0) / Real_t(3.0) ;
}
#if USE_MPI
MPI_Allreduce(&gnewdt, &newdt, 1,
((sizeof(Real_t) == 4) ? MPI_FLOAT : MPI_DOUBLE),
MPI_MIN, MPI_COMM_WORLD) ;
#else
newdt = gnewdt;
#endif
ratio = newdt / olddt ;
if (ratio >= Real_t(1.0)) {
if (ratio < domain.deltatimemultlb()) {
newdt = olddt ;
}
else if (ratio > domain.deltatimemultub()) {
newdt = olddt*domain.deltatimemultub() ;
}
}
if (newdt > domain.dtmax()) {
newdt = domain.dtmax() ;
}
domain.deltatime() = newdt ;
}
/* TRY TO PREVENT VERY SMALL SCALING ON THE NEXT CYCLE */
if ((targetdt > domain.deltatime()) &&
(targetdt < (Real_t(4.0) * domain.deltatime() / Real_t(3.0))) ) {
targetdt = Real_t(2.0) * domain.deltatime() / Real_t(3.0) ;
}
if (targetdt < domain.deltatime()) {
domain.deltatime() = targetdt ;
}
domain.time() += domain.deltatime() ;
++domain.cycle() ;
}
/******************************************/
static inline
void CollectDomainNodesToElemNodes(Domain &domain,
const Index_t* elemToNode,
Real_t elemX[8],
Real_t elemY[8],
Real_t elemZ[8])
{
Index_t nd0i = elemToNode[0] ;
Index_t nd1i = elemToNode[1] ;
Index_t nd2i = elemToNode[2] ;
Index_t nd3i = elemToNode[3] ;
Index_t nd4i = elemToNode[4] ;
Index_t nd5i = elemToNode[5] ;
Index_t nd6i = elemToNode[6] ;
Index_t nd7i = elemToNode[7] ;
elemX[0] = domain.x(nd0i);
elemX[1] = domain.x(nd1i);
elemX[2] = domain.x(nd2i);
elemX[3] = domain.x(nd3i);
elemX[4] = domain.x(nd4i);
elemX[5] = domain.x(nd5i);
elemX[6] = domain.x(nd6i);
elemX[7] = domain.x(nd7i);
elemY[0] = domain.y(nd0i);
elemY[1] = domain.y(nd1i);
elemY[2] = domain.y(nd2i);
elemY[3] = domain.y(nd3i);
elemY[4] = domain.y(nd4i);
elemY[5] = domain.y(nd5i);
elemY[6] = domain.y(nd6i);
elemY[7] = domain.y(nd7i);
elemZ[0] = domain.z(nd0i);
elemZ[1] = domain.z(nd1i);
elemZ[2] = domain.z(nd2i);
elemZ[3] = domain.z(nd3i);
elemZ[4] = domain.z(nd4i);
elemZ[5] = domain.z(nd5i);
elemZ[6] = domain.z(nd6i);
elemZ[7] = domain.z(nd7i);
}
/******************************************/
static inline
void InitStressTermsForElems(Domain &domain,
Real_t *sigxx, Real_t *sigyy, Real_t *sigzz,
Index_t numElem)
{
//
// pull in the stresses appropriate to the hydro integration
//
#pragma omp parallel for firstprivate(numElem)
for (Index_t i = 0 ; i < numElem ; ++i){
sigxx[i] = sigyy[i] = sigzz[i] = - domain.p(i) - domain.q(i) ;
}
}
/******************************************/
static inline
void CalcElemShapeFunctionDerivatives( Real_t const x[],
Real_t const y[],
Real_t const z[],
Real_t b[][8],
Real_t* const volume )
{
const Real_t x0 = x[0] ; const Real_t x1 = x[1] ;
const Real_t x2 = x[2] ; const Real_t x3 = x[3] ;
const Real_t x4 = x[4] ; const Real_t x5 = x[5] ;
const Real_t x6 = x[6] ; const Real_t x7 = x[7] ;
const Real_t y0 = y[0] ; const Real_t y1 = y[1] ;
const Real_t y2 = y[2] ; const Real_t y3 = y[3] ;
const Real_t y4 = y[4] ; const Real_t y5 = y[5] ;
const Real_t y6 = y[6] ; const Real_t y7 = y[7] ;
const Real_t z0 = z[0] ; const Real_t z1 = z[1] ;
const Real_t z2 = z[2] ; const Real_t z3 = z[3] ;
const Real_t z4 = z[4] ; const Real_t z5 = z[5] ;
const Real_t z6 = z[6] ; const Real_t z7 = z[7] ;
Real_t fjxxi, fjxet, fjxze;
Real_t fjyxi, fjyet, fjyze;
Real_t fjzxi, fjzet, fjzze;
Real_t cjxxi, cjxet, cjxze;
Real_t cjyxi, cjyet, cjyze;
Real_t cjzxi, cjzet, cjzze;
fjxxi = Real_t(.125) * ( (x6-x0) + (x5-x3) - (x7-x1) - (x4-x2) );
fjxet = Real_t(.125) * ( (x6-x0) - (x5-x3) + (x7-x1) - (x4-x2) );
fjxze = Real_t(.125) * ( (x6-x0) + (x5-x3) + (x7-x1) + (x4-x2) );
fjyxi = Real_t(.125) * ( (y6-y0) + (y5-y3) - (y7-y1) - (y4-y2) );
fjyet = Real_t(.125) * ( (y6-y0) - (y5-y3) + (y7-y1) - (y4-y2) );
fjyze = Real_t(.125) * ( (y6-y0) + (y5-y3) + (y7-y1) + (y4-y2) );
fjzxi = Real_t(.125) * ( (z6-z0) + (z5-z3) - (z7-z1) - (z4-z2) );
fjzet = Real_t(.125) * ( (z6-z0) - (z5-z3) + (z7-z1) - (z4-z2) );
fjzze = Real_t(.125) * ( (z6-z0) + (z5-z3) + (z7-z1) + (z4-z2) );
/* compute cofactors */
cjxxi = (fjyet * fjzze) - (fjzet * fjyze);
cjxet = - (fjyxi * fjzze) + (fjzxi * fjyze);
cjxze = (fjyxi * fjzet) - (fjzxi * fjyet);
cjyxi = - (fjxet * fjzze) + (fjzet * fjxze);
cjyet = (fjxxi * fjzze) - (fjzxi * fjxze);
cjyze = - (fjxxi * fjzet) + (fjzxi * fjxet);
cjzxi = (fjxet * fjyze) - (fjyet * fjxze);
cjzet = - (fjxxi * fjyze) + (fjyxi * fjxze);
cjzze = (fjxxi * fjyet) - (fjyxi * fjxet);
/* calculate partials :
this need only be done for l = 0,1,2,3 since , by symmetry ,
(6,7,4,5) = - (0,1,2,3) .
*/
b[0][0] = - cjxxi - cjxet - cjxze;
b[0][1] = cjxxi - cjxet - cjxze;
b[0][2] = cjxxi + cjxet - cjxze;
b[0][3] = - cjxxi + cjxet - cjxze;
b[0][4] = -b[0][2];
b[0][5] = -b[0][3];
b[0][6] = -b[0][0];
b[0][7] = -b[0][1];
b[1][0] = - cjyxi - cjyet - cjyze;
b[1][1] = cjyxi - cjyet - cjyze;
b[1][2] = cjyxi + cjyet - cjyze;
b[1][3] = - cjyxi + cjyet - cjyze;
b[1][4] = -b[1][2];
b[1][5] = -b[1][3];
b[1][6] = -b[1][0];
b[1][7] = -b[1][1];
b[2][0] = - cjzxi - cjzet - cjzze;
b[2][1] = cjzxi - cjzet - cjzze;
b[2][2] = cjzxi + cjzet - cjzze;
b[2][3] = - cjzxi + cjzet - cjzze;
b[2][4] = -b[2][2];
b[2][5] = -b[2][3];
b[2][6] = -b[2][0];
b[2][7] = -b[2][1];
/* calculate jacobian determinant (volume) */
*volume = Real_t(8.) * ( fjxet * cjxet + fjyet * cjyet + fjzet * cjzet);
}
/******************************************/
static inline
void SumElemFaceNormal(Real_t *normalX0, Real_t *normalY0, Real_t *normalZ0,
Real_t *normalX1, Real_t *normalY1, Real_t *normalZ1,
Real_t *normalX2, Real_t *normalY2, Real_t *normalZ2,
Real_t *normalX3, Real_t *normalY3, Real_t *normalZ3,
const Real_t x0, const Real_t y0, const Real_t z0,
const Real_t x1, const Real_t y1, const Real_t z1,
const Real_t x2, const Real_t y2, const Real_t z2,
const Real_t x3, const Real_t y3, const Real_t z3)
{
Real_t bisectX0 = Real_t(0.5) * (x3 + x2 - x1 - x0);
Real_t bisectY0 = Real_t(0.5) * (y3 + y2 - y1 - y0);
Real_t bisectZ0 = Real_t(0.5) * (z3 + z2 - z1 - z0);
Real_t bisectX1 = Real_t(0.5) * (x2 + x1 - x3 - x0);
Real_t bisectY1 = Real_t(0.5) * (y2 + y1 - y3 - y0);
Real_t bisectZ1 = Real_t(0.5) * (z2 + z1 - z3 - z0);
Real_t areaX = Real_t(0.25) * (bisectY0 * bisectZ1 - bisectZ0 * bisectY1);
Real_t areaY = Real_t(0.25) * (bisectZ0 * bisectX1 - bisectX0 * bisectZ1);
Real_t areaZ = Real_t(0.25) * (bisectX0 * bisectY1 - bisectY0 * bisectX1);
*normalX0 += areaX;
*normalX1 += areaX;
*normalX2 += areaX;
*normalX3 += areaX;
*normalY0 += areaY;
*normalY1 += areaY;
*normalY2 += areaY;
*normalY3 += areaY;
*normalZ0 += areaZ;
*normalZ1 += areaZ;
*normalZ2 += areaZ;
*normalZ3 += areaZ;
}
/******************************************/
static inline
void CalcElemNodeNormals(Real_t pfx[8],
Real_t pfy[8],
Real_t pfz[8],
const Real_t x[8],
const Real_t y[8],
const Real_t z[8])
{
for (Index_t i = 0 ; i < 8 ; ++i) {
pfx[i] = Real_t(0.0);
pfy[i] = Real_t(0.0);
pfz[i] = Real_t(0.0);
}
/* evaluate face one: nodes 0, 1, 2, 3 */
SumElemFaceNormal(&pfx[0], &pfy[0], &pfz[0],
&pfx[1], &pfy[1], &pfz[1],
&pfx[2], &pfy[2], &pfz[2],
&pfx[3], &pfy[3], &pfz[3],
x[0], y[0], z[0], x[1], y[1], z[1],
x[2], y[2], z[2], x[3], y[3], z[3]);
/* evaluate face two: nodes 0, 4, 5, 1 */
SumElemFaceNormal(&pfx[0], &pfy[0], &pfz[0],
&pfx[4], &pfy[4], &pfz[4],
&pfx[5], &pfy[5], &pfz[5],
&pfx[1], &pfy[1], &pfz[1],
x[0], y[0], z[0], x[4], y[4], z[4],
x[5], y[5], z[5], x[1], y[1], z[1]);
/* evaluate face three: nodes 1, 5, 6, 2 */
SumElemFaceNormal(&pfx[1], &pfy[1], &pfz[1],
&pfx[5], &pfy[5], &pfz[5],
&pfx[6], &pfy[6], &pfz[6],
&pfx[2], &pfy[2], &pfz[2],
x[1], y[1], z[1], x[5], y[5], z[5],
x[6], y[6], z[6], x[2], y[2], z[2]);
/* evaluate face four: nodes 2, 6, 7, 3 */
SumElemFaceNormal(&pfx[2], &pfy[2], &pfz[2],
&pfx[6], &pfy[6], &pfz[6],
&pfx[7], &pfy[7], &pfz[7],
&pfx[3], &pfy[3], &pfz[3],
x[2], y[2], z[2], x[6], y[6], z[6],
x[7], y[7], z[7], x[3], y[3], z[3]);
/* evaluate face five: nodes 3, 7, 4, 0 */
SumElemFaceNormal(&pfx[3], &pfy[3], &pfz[3],
&pfx[7], &pfy[7], &pfz[7],
&pfx[4], &pfy[4], &pfz[4],
&pfx[0], &pfy[0], &pfz[0],
x[3], y[3], z[3], x[7], y[7], z[7],
x[4], y[4], z[4], x[0], y[0], z[0]);
/* evaluate face six: nodes 4, 7, 6, 5 */
SumElemFaceNormal(&pfx[4], &pfy[4], &pfz[4],
&pfx[7], &pfy[7], &pfz[7],
&pfx[6], &pfy[6], &pfz[6],
&pfx[5], &pfy[5], &pfz[5],
x[4], y[4], z[4], x[7], y[7], z[7],
x[6], y[6], z[6], x[5], y[5], z[5]);
}
/******************************************/
static inline
void SumElemStressesToNodeForces( const Real_t B[][8],
const Real_t stress_xx,
const Real_t stress_yy,
const Real_t stress_zz,
Real_t fx[], Real_t fy[], Real_t fz[] )
{
for(Index_t i = 0; i < 8; i++) {
fx[i] = -( stress_xx * B[0][i] );
fy[i] = -( stress_yy * B[1][i] );
fz[i] = -( stress_zz * B[2][i] );
}
}
/******************************************/
static inline
void IntegrateStressForElems( Domain &domain,
Real_t *sigxx, Real_t *sigyy, Real_t *sigzz,
Real_t *determ, Index_t numElem, Index_t numNode)
{
#if _OPENMP
Index_t numthreads = omp_get_max_threads();
#else
Index_t numthreads = 1;
#endif
Index_t numElem8 = numElem * 8 ;
Real_t *fx_elem;
Real_t *fy_elem;
Real_t *fz_elem;
Real_t fx_local[8] ;
Real_t fy_local[8] ;
Real_t fz_local[8] ;
if (numthreads > 1) {
fx_elem = Allocate<Real_t>(numElem8) ;
fy_elem = Allocate<Real_t>(numElem8) ;
fz_elem = Allocate<Real_t>(numElem8) ;
}
// loop over all elements
#pragma omp parallel for firstprivate(numElem)
for( Index_t k=0 ; k<numElem ; ++k )
{
const Index_t* const elemToNode = domain.nodelist(k);
Real_t B[3][8] ;// shape function derivatives
Real_t x_local[8] ;
Real_t y_local[8] ;
Real_t z_local[8] ;
// get nodal coordinates from global arrays and copy into local arrays.
CollectDomainNodesToElemNodes(domain, elemToNode, x_local, y_local, z_local);
// Volume calculation involves extra work for numerical consistency
CalcElemShapeFunctionDerivatives(x_local, y_local, z_local,
B, &determ[k]);
CalcElemNodeNormals( B[0] , B[1], B[2],
x_local, y_local, z_local );
if (numthreads > 1) {
// Eliminate thread writing conflicts at the nodes by giving
// each element its own copy to write to
SumElemStressesToNodeForces( B, sigxx[k], sigyy[k], sigzz[k],
&fx_elem[k*8],
&fy_elem[k*8],
&fz_elem[k*8] ) ;
}
else {
SumElemStressesToNodeForces( B, sigxx[k], sigyy[k], sigzz[k],
fx_local, fy_local, fz_local ) ;
// copy nodal force contributions to global force arrray.
for( Index_t lnode=0 ; lnode<8 ; ++lnode ) {
Index_t gnode = elemToNode[lnode];
domain.fx(gnode) += fx_local[lnode];
domain.fy(gnode) += fy_local[lnode];
domain.fz(gnode) += fz_local[lnode];
}
}
}
if (numthreads > 1) {
// If threaded, then we need to copy the data out of the temporary
// arrays used above into the final forces field
#pragma omp parallel for firstprivate(numNode)
for( Index_t gnode=0 ; gnode<numNode ; ++gnode )
{
Index_t count = domain.nodeElemCount(gnode) ;
Index_t *cornerList = domain.nodeElemCornerList(gnode) ;
Real_t fx_tmp = Real_t(0.0) ;
Real_t fy_tmp = Real_t(0.0) ;
Real_t fz_tmp = Real_t(0.0) ;
for (Index_t i=0 ; i < count ; ++i) {
Index_t ielem = cornerList[i] ;
fx_tmp += fx_elem[ielem] ;
fy_tmp += fy_elem[ielem] ;
fz_tmp += fz_elem[ielem] ;
}
domain.fx(gnode) = fx_tmp ;
domain.fy(gnode) = fy_tmp ;
domain.fz(gnode) = fz_tmp ;
}
Release(&fz_elem) ;
Release(&fy_elem) ;
Release(&fx_elem) ;
}
}
/******************************************/
static inline
void VoluDer(const Real_t x0, const Real_t x1, const Real_t x2,
const Real_t x3, const Real_t x4, const Real_t x5,
const Real_t y0, const Real_t y1, const Real_t y2,
const Real_t y3, const Real_t y4, const Real_t y5,
const Real_t z0, const Real_t z1, const Real_t z2,
const Real_t z3, const Real_t z4, const Real_t z5,
Real_t* dvdx, Real_t* dvdy, Real_t* dvdz)
{
const Real_t twelfth = Real_t(1.0) / Real_t(12.0) ;
*dvdx =
(y1 + y2) * (z0 + z1) - (y0 + y1) * (z1 + z2) +
(y0 + y4) * (z3 + z4) - (y3 + y4) * (z0 + z4) -
(y2 + y5) * (z3 + z5) + (y3 + y5) * (z2 + z5);
*dvdy =
- (x1 + x2) * (z0 + z1) + (x0 + x1) * (z1 + z2) -
(x0 + x4) * (z3 + z4) + (x3 + x4) * (z0 + z4) +
(x2 + x5) * (z3 + z5) - (x3 + x5) * (z2 + z5);
*dvdz =
- (y1 + y2) * (x0 + x1) + (y0 + y1) * (x1 + x2) -
(y0 + y4) * (x3 + x4) + (y3 + y4) * (x0 + x4) +
(y2 + y5) * (x3 + x5) - (y3 + y5) * (x2 + x5);
*dvdx *= twelfth;
*dvdy *= twelfth;
*dvdz *= twelfth;
}
/******************************************/
static inline
void CalcElemVolumeDerivative(Real_t dvdx[8],
Real_t dvdy[8],
Real_t dvdz[8],
const Real_t x[8],
const Real_t y[8],
const Real_t z[8])
{
VoluDer(x[1], x[2], x[3], x[4], x[5], x[7],
y[1], y[2], y[3], y[4], y[5], y[7],
z[1], z[2], z[3], z[4], z[5], z[7],
&dvdx[0], &dvdy[0], &dvdz[0]);
VoluDer(x[0], x[1], x[2], x[7], x[4], x[6],
y[0], y[1], y[2], y[7], y[4], y[6],
z[0], z[1], z[2], z[7], z[4], z[6],
&dvdx[3], &dvdy[3], &dvdz[3]);
VoluDer(x[3], x[0], x[1], x[6], x[7], x[5],
y[3], y[0], y[1], y[6], y[7], y[5],
z[3], z[0], z[1], z[6], z[7], z[5],
&dvdx[2], &dvdy[2], &dvdz[2]);
VoluDer(x[2], x[3], x[0], x[5], x[6], x[4],
y[2], y[3], y[0], y[5], y[6], y[4],
z[2], z[3], z[0], z[5], z[6], z[4],
&dvdx[1], &dvdy[1], &dvdz[1]);
VoluDer(x[7], x[6], x[5], x[0], x[3], x[1],
y[7], y[6], y[5], y[0], y[3], y[1],
z[7], z[6], z[5], z[0], z[3], z[1],
&dvdx[4], &dvdy[4], &dvdz[4]);
VoluDer(x[4], x[7], x[6], x[1], x[0], x[2],
y[4], y[7], y[6], y[1], y[0], y[2],
z[4], z[7], z[6], z[1], z[0], z[2],
&dvdx[5], &dvdy[5], &dvdz[5]);
VoluDer(x[5], x[4], x[7], x[2], x[1], x[3],
y[5], y[4], y[7], y[2], y[1], y[3],
z[5], z[4], z[7], z[2], z[1], z[3],
&dvdx[6], &dvdy[6], &dvdz[6]);
VoluDer(x[6], x[5], x[4], x[3], x[2], x[0],
y[6], y[5], y[4], y[3], y[2], y[0],
z[6], z[5], z[4], z[3], z[2], z[0],
&dvdx[7], &dvdy[7], &dvdz[7]);
}
/******************************************/
static inline
void CalcElemFBHourglassForce(Real_t *xd, Real_t *yd, Real_t *zd, Real_t hourgam[][4],
Real_t coefficient,
Real_t *hgfx, Real_t *hgfy, Real_t *hgfz )
{
Real_t hxx[4];
for(Index_t i = 0; i < 4; i++) {
hxx[i] = hourgam[0][i] * xd[0] + hourgam[1][i] * xd[1] +
hourgam[2][i] * xd[2] + hourgam[3][i] * xd[3] +
hourgam[4][i] * xd[4] + hourgam[5][i] * xd[5] +
hourgam[6][i] * xd[6] + hourgam[7][i] * xd[7];
}
for(Index_t i = 0; i < 8; i++) {
hgfx[i] = coefficient *
(hourgam[i][0] * hxx[0] + hourgam[i][1] * hxx[1] +
hourgam[i][2] * hxx[2] + hourgam[i][3] * hxx[3]);
}
for(Index_t i = 0; i < 4; i++) {
hxx[i] = hourgam[0][i] * yd[0] + hourgam[1][i] * yd[1] +
hourgam[2][i] * yd[2] + hourgam[3][i] * yd[3] +
hourgam[4][i] * yd[4] + hourgam[5][i] * yd[5] +
hourgam[6][i] * yd[6] + hourgam[7][i] * yd[7];
}
for(Index_t i = 0; i < 8; i++) {
hgfy[i] = coefficient *
(hourgam[i][0] * hxx[0] + hourgam[i][1] * hxx[1] +
hourgam[i][2] * hxx[2] + hourgam[i][3] * hxx[3]);
}
for(Index_t i = 0; i < 4; i++) {
hxx[i] = hourgam[0][i] * zd[0] + hourgam[1][i] * zd[1] +
hourgam[2][i] * zd[2] + hourgam[3][i] * zd[3] +
hourgam[4][i] * zd[4] + hourgam[5][i] * zd[5] +
hourgam[6][i] * zd[6] + hourgam[7][i] * zd[7];
}
for(Index_t i = 0; i < 8; i++) {
hgfz[i] = coefficient *
(hourgam[i][0] * hxx[0] + hourgam[i][1] * hxx[1] +
hourgam[i][2] * hxx[2] + hourgam[i][3] * hxx[3]);
}
}
/******************************************/
static inline
void CalcFBHourglassForceForElems( Domain &domain,
Real_t *determ,
Real_t *x8n, Real_t *y8n, Real_t *z8n,
Real_t *dvdx, Real_t *dvdy, Real_t *dvdz,
Real_t hourg, Index_t numElem,
Index_t numNode)
{
#if _OPENMP
Index_t numthreads = omp_get_max_threads();
#else
Index_t numthreads = 1;
#endif
/*************************************************
*
* FUNCTION: Calculates the Flanagan-Belytschko anti-hourglass
* force.
*
*************************************************/
Index_t numElem8 = numElem * 8 ;
Real_t *fx_elem;
Real_t *fy_elem;
Real_t *fz_elem;
if(numthreads > 1) {
fx_elem = Allocate<Real_t>(numElem8) ;
fy_elem = Allocate<Real_t>(numElem8) ;
fz_elem = Allocate<Real_t>(numElem8) ;
}
Real_t gamma[4][8];
gamma[0][0] = Real_t( 1.);
gamma[0][1] = Real_t( 1.);
gamma[0][2] = Real_t(-1.);
gamma[0][3] = Real_t(-1.);
gamma[0][4] = Real_t(-1.);
gamma[0][5] = Real_t(-1.);
gamma[0][6] = Real_t( 1.);
gamma[0][7] = Real_t( 1.);
gamma[1][0] = Real_t( 1.);
gamma[1][1] = Real_t(-1.);
gamma[1][2] = Real_t(-1.);
gamma[1][3] = Real_t( 1.);
gamma[1][4] = Real_t(-1.);
gamma[1][5] = Real_t( 1.);
gamma[1][6] = Real_t( 1.);
gamma[1][7] = Real_t(-1.);
gamma[2][0] = Real_t( 1.);
gamma[2][1] = Real_t(-1.);
gamma[2][2] = Real_t( 1.);
gamma[2][3] = Real_t(-1.);
gamma[2][4] = Real_t( 1.);
gamma[2][5] = Real_t(-1.);
gamma[2][6] = Real_t( 1.);
gamma[2][7] = Real_t(-1.);
gamma[3][0] = Real_t(-1.);
gamma[3][1] = Real_t( 1.);
gamma[3][2] = Real_t(-1.);
gamma[3][3] = Real_t( 1.);
gamma[3][4] = Real_t( 1.);
gamma[3][5] = Real_t(-1.);
gamma[3][6] = Real_t( 1.);
gamma[3][7] = Real_t(-1.);
/*************************************************/
/* compute the hourglass modes */
#pragma omp parallel for firstprivate(numElem, hourg)
for(Index_t i2=0;i2<numElem;++i2){
Real_t *fx_local, *fy_local, *fz_local ;
Real_t hgfx[8], hgfy[8], hgfz[8] ;
Real_t coefficient;
Real_t hourgam[8][4];
Real_t xd1[8], yd1[8], zd1[8] ;
const Index_t *elemToNode = domain.nodelist(i2);
Index_t i3=8*i2;
Real_t volinv=Real_t(1.0)/determ[i2];
Real_t ss1, mass1, volume13 ;
for(Index_t i1=0;i1<4;++i1){
Real_t hourmodx =
x8n[i3] * gamma[i1][0] + x8n[i3+1] * gamma[i1][1] +
x8n[i3+2] * gamma[i1][2] + x8n[i3+3] * gamma[i1][3] +
x8n[i3+4] * gamma[i1][4] + x8n[i3+5] * gamma[i1][5] +
x8n[i3+6] * gamma[i1][6] + x8n[i3+7] * gamma[i1][7];
Real_t hourmody =
y8n[i3] * gamma[i1][0] + y8n[i3+1] * gamma[i1][1] +
y8n[i3+2] * gamma[i1][2] + y8n[i3+3] * gamma[i1][3] +
y8n[i3+4] * gamma[i1][4] + y8n[i3+5] * gamma[i1][5] +
y8n[i3+6] * gamma[i1][6] + y8n[i3+7] * gamma[i1][7];
Real_t hourmodz =
z8n[i3] * gamma[i1][0] + z8n[i3+1] * gamma[i1][1] +
z8n[i3+2] * gamma[i1][2] + z8n[i3+3] * gamma[i1][3] +
z8n[i3+4] * gamma[i1][4] + z8n[i3+5] * gamma[i1][5] +
z8n[i3+6] * gamma[i1][6] + z8n[i3+7] * gamma[i1][7];
hourgam[0][i1] = gamma[i1][0] - volinv*(dvdx[i3 ] * hourmodx +
dvdy[i3 ] * hourmody +
dvdz[i3 ] * hourmodz );
hourgam[1][i1] = gamma[i1][1] - volinv*(dvdx[i3+1] * hourmodx +
dvdy[i3+1] * hourmody +
dvdz[i3+1] * hourmodz );
hourgam[2][i1] = gamma[i1][2] - volinv*(dvdx[i3+2] * hourmodx +
dvdy[i3+2] * hourmody +
dvdz[i3+2] * hourmodz );
hourgam[3][i1] = gamma[i1][3] - volinv*(dvdx[i3+3] * hourmodx +
dvdy[i3+3] * hourmody +
dvdz[i3+3] * hourmodz );
hourgam[4][i1] = gamma[i1][4] - volinv*(dvdx[i3+4] * hourmodx +
dvdy[i3+4] * hourmody +
dvdz[i3+4] * hourmodz );
hourgam[5][i1] = gamma[i1][5] - volinv*(dvdx[i3+5] * hourmodx +
dvdy[i3+5] * hourmody +
dvdz[i3+5] * hourmodz );
hourgam[6][i1] = gamma[i1][6] - volinv*(dvdx[i3+6] * hourmodx +
dvdy[i3+6] * hourmody +
dvdz[i3+6] * hourmodz );
hourgam[7][i1] = gamma[i1][7] - volinv*(dvdx[i3+7] * hourmodx +
dvdy[i3+7] * hourmody +
dvdz[i3+7] * hourmodz );
}
/* compute forces */
/* store forces into h arrays (force arrays) */
ss1=domain.ss(i2);
mass1=domain.elemMass(i2);
volume13=CBRT(determ[i2]);
Index_t n0si2 = elemToNode[0];
Index_t n1si2 = elemToNode[1];
Index_t n2si2 = elemToNode[2];
Index_t n3si2 = elemToNode[3];
Index_t n4si2 = elemToNode[4];
Index_t n5si2 = elemToNode[5];
Index_t n6si2 = elemToNode[6];
Index_t n7si2 = elemToNode[7];
xd1[0] = domain.xd(n0si2);
xd1[1] = domain.xd(n1si2);
xd1[2] = domain.xd(n2si2);
xd1[3] = domain.xd(n3si2);
xd1[4] = domain.xd(n4si2);
xd1[5] = domain.xd(n5si2);
xd1[6] = domain.xd(n6si2);
xd1[7] = domain.xd(n7si2);
yd1[0] = domain.yd(n0si2);
yd1[1] = domain.yd(n1si2);
yd1[2] = domain.yd(n2si2);
yd1[3] = domain.yd(n3si2);
yd1[4] = domain.yd(n4si2);
yd1[5] = domain.yd(n5si2);
yd1[6] = domain.yd(n6si2);
yd1[7] = domain.yd(n7si2);
zd1[0] = domain.zd(n0si2);
zd1[1] = domain.zd(n1si2);
zd1[2] = domain.zd(n2si2);
zd1[3] = domain.zd(n3si2);
zd1[4] = domain.zd(n4si2);
zd1[5] = domain.zd(n5si2);
zd1[6] = domain.zd(n6si2);
zd1[7] = domain.zd(n7si2);
coefficient = - hourg * Real_t(0.01) * ss1 * mass1 / volume13;
CalcElemFBHourglassForce(xd1,yd1,zd1,
hourgam,
coefficient, hgfx, hgfy, hgfz);
// With the threaded version, we write into local arrays per elem
// so we don't have to worry about race conditions
if (numthreads > 1) {
fx_local = &fx_elem[i3] ;
fx_local[0] = hgfx[0];
fx_local[1] = hgfx[1];
fx_local[2] = hgfx[2];
fx_local[3] = hgfx[3];
fx_local[4] = hgfx[4];
fx_local[5] = hgfx[5];
fx_local[6] = hgfx[6];
fx_local[7] = hgfx[7];
fy_local = &fy_elem[i3] ;
fy_local[0] = hgfy[0];
fy_local[1] = hgfy[1];
fy_local[2] = hgfy[2];
fy_local[3] = hgfy[3];
fy_local[4] = hgfy[4];
fy_local[5] = hgfy[5];
fy_local[6] = hgfy[6];
fy_local[7] = hgfy[7];
fz_local = &fz_elem[i3] ;
fz_local[0] = hgfz[0];
fz_local[1] = hgfz[1];
fz_local[2] = hgfz[2];
fz_local[3] = hgfz[3];
fz_local[4] = hgfz[4];
fz_local[5] = hgfz[5];
fz_local[6] = hgfz[6];
fz_local[7] = hgfz[7];
}
else {
domain.fx(n0si2) += hgfx[0];
domain.fy(n0si2) += hgfy[0];
domain.fz(n0si2) += hgfz[0];
domain.fx(n1si2) += hgfx[1];
domain.fy(n1si2) += hgfy[1];
domain.fz(n1si2) += hgfz[1];
domain.fx(n2si2) += hgfx[2];
domain.fy(n2si2) += hgfy[2];
domain.fz(n2si2) += hgfz[2];
domain.fx(n3si2) += hgfx[3];
domain.fy(n3si2) += hgfy[3];
domain.fz(n3si2) += hgfz[3];
domain.fx(n4si2) += hgfx[4];
domain.fy(n4si2) += hgfy[4];
domain.fz(n4si2) += hgfz[4];
domain.fx(n5si2) += hgfx[5];
domain.fy(n5si2) += hgfy[5];
domain.fz(n5si2) += hgfz[5];
domain.fx(n6si2) += hgfx[6];
domain.fy(n6si2) += hgfy[6];
domain.fz(n6si2) += hgfz[6];
domain.fx(n7si2) += hgfx[7];
domain.fy(n7si2) += hgfy[7];
domain.fz(n7si2) += hgfz[7];
}
}
if (numthreads > 1) {
// Collect the data from the local arrays into the final force arrays
#pragma omp parallel for firstprivate(numNode)
for( Index_t gnode=0 ; gnode<numNode ; ++gnode )
{
Index_t count = domain.nodeElemCount(gnode) ;
Index_t *cornerList = domain.nodeElemCornerList(gnode) ;
Real_t fx_tmp = Real_t(0.0) ;
Real_t fy_tmp = Real_t(0.0) ;
Real_t fz_tmp = Real_t(0.0) ;
for (Index_t i=0 ; i < count ; ++i) {
Index_t ielem = cornerList[i] ;
fx_tmp += fx_elem[ielem] ;
fy_tmp += fy_elem[ielem] ;
fz_tmp += fz_elem[ielem] ;
}
domain.fx(gnode) += fx_tmp ;
domain.fy(gnode) += fy_tmp ;
domain.fz(gnode) += fz_tmp ;
}
Release(&fz_elem) ;
Release(&fy_elem) ;
Release(&fx_elem) ;
}
}
/******************************************/
static inline
void CalcHourglassControlForElems(Domain& domain,
Real_t determ[], Real_t hgcoef)
{
Index_t numElem = domain.numElem() ;
Index_t numElem8 = numElem * 8 ;