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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include "modify.h"
#include <cstring>
#include "style_compute.h"
#include "style_fix.h"
#include "atom.h"
#include "comm.h"
#include "fix.h"
#include "compute.h"
#include "group.h"
#include "update.h"
#include "domain.h"
#include "region.h"
#include "input.h"
#include "variable.h"
#include "memory.h"
#include "error.h"
#include "utils.h"
using namespace LAMMPS_NS;
using namespace FixConst;
#define DELTA 4
#define BIG 1.0e20
#define NEXCEPT 7 // change when add to exceptions in add_fix()
/* ---------------------------------------------------------------------- */
Modify::Modify(LAMMPS *lmp) : Pointers(lmp)
{
nfix = maxfix = 0;
n_initial_integrate = n_post_integrate = 0;
n_pre_exchange = n_pre_neighbor = n_post_neighbor = 0;
n_pre_force = n_pre_reverse = n_post_force = 0;
n_final_integrate = n_end_of_step = n_thermo_energy = 0;
n_thermo_energy_atom = 0;
n_initial_integrate_respa = n_post_integrate_respa = 0;
n_pre_force_respa = n_post_force_respa = n_final_integrate_respa = 0;
n_min_pre_exchange = n_min_pre_force = n_min_pre_reverse = 0;
n_min_post_force = n_min_energy = 0;
fix = NULL;
fmask = NULL;
list_initial_integrate = list_post_integrate = NULL;
list_pre_exchange = list_pre_neighbor = list_post_neighbor = NULL;
list_pre_force = list_pre_reverse = list_post_force = NULL;
list_final_integrate = list_end_of_step = NULL;
list_thermo_energy = list_thermo_energy_atom = NULL;
list_initial_integrate_respa = list_post_integrate_respa = NULL;
list_pre_force_respa = list_post_force_respa = NULL;
list_final_integrate_respa = NULL;
list_min_pre_exchange = list_min_pre_neighbor = list_min_post_neighbor = NULL;
list_min_pre_force = list_min_pre_reverse = list_min_post_force = NULL;
list_min_energy = NULL;
end_of_step_every = NULL;
list_timeflag = NULL;
nfix_restart_global = 0;
id_restart_global = style_restart_global = NULL;
state_restart_global = NULL;
used_restart_global = NULL;
nfix_restart_peratom = 0;
id_restart_peratom = style_restart_peratom = NULL;
index_restart_peratom = used_restart_peratom = NULL;
ncompute = maxcompute = 0;
compute = NULL;
create_factories();
}
void _noopt Modify::create_factories()
{
// fill map with fixes listed in style_fix.h
fix_map = new FixCreatorMap();
#define FIX_CLASS
#define FixStyle(key,Class) \
(*fix_map)[#key] = &fix_creator<Class>;
#include "style_fix.h"
#undef FixStyle
#undef FIX_CLASS
// fill map with computes listed in style_compute.h
compute_map = new ComputeCreatorMap();
#define COMPUTE_CLASS
#define ComputeStyle(key,Class) \
(*compute_map)[#key] = &compute_creator<Class>;
#include "style_compute.h"
#undef ComputeStyle
#undef COMPUTE_CLASS
}
/* ---------------------------------------------------------------------- */
Modify::~Modify()
{
// delete all fixes
// do it via delete_fix() so callbacks in Atom are also updated correctly
while (nfix) delete_fix(0);
memory->sfree(fix);
memory->destroy(fmask);
// delete all computes
for (int i = 0; i < ncompute; i++) delete compute[i];
memory->sfree(compute);
delete [] list_initial_integrate;
delete [] list_post_integrate;
delete [] list_pre_exchange;
delete [] list_pre_neighbor;
delete [] list_post_neighbor;
delete [] list_pre_force;
delete [] list_pre_reverse;
delete [] list_post_force;
delete [] list_final_integrate;
delete [] list_end_of_step;
delete [] list_thermo_energy;
delete [] list_thermo_energy_atom;
delete [] list_initial_integrate_respa;
delete [] list_post_integrate_respa;
delete [] list_pre_force_respa;
delete [] list_post_force_respa;
delete [] list_final_integrate_respa;
delete [] list_min_pre_exchange;
delete [] list_min_pre_neighbor;
delete [] list_min_post_neighbor;
delete [] list_min_pre_force;
delete [] list_min_pre_reverse;
delete [] list_min_post_force;
delete [] list_min_energy;
delete [] end_of_step_every;
delete [] list_timeflag;
restart_deallocate(0);
delete compute_map;
delete fix_map;
}
/* ----------------------------------------------------------------------
initialize all fixes and computes
------------------------------------------------------------------------- */
void Modify::init()
{
int i,j;
// delete storage of restart info since it is not valid after 1st run
restart_deallocate(1);
// init each compute
// set invoked_scalar,vector,etc to -1 to force new run to re-compute them
// add initial timestep to all computes that store invocation times
// since any of them may be invoked by initial thermo
// do not clear out invocation times stored within a compute,
// b/c some may be holdovers from previous run, like for ave fixes
for (i = 0; i < ncompute; i++) {
compute[i]->init();
compute[i]->invoked_scalar = -1;
compute[i]->invoked_vector = -1;
compute[i]->invoked_array = -1;
compute[i]->invoked_peratom = -1;
compute[i]->invoked_local = -1;
}
addstep_compute_all(update->ntimestep);
// init each fix
// should not need to come before compute init
// used to b/c temperature computes called fix->dof() in their init,
// and fix rigid required its own init before its dof() could be called,
// but computes now do their DOF in setup()
for (i = 0; i < nfix; i++) fix[i]->init();
// set global flag if any fix has its restart_pbc flag set
restart_pbc_any = 0;
for (i = 0; i < nfix; i++)
if (fix[i]->restart_pbc) restart_pbc_any = 1;
// create lists of fixes to call at each stage of run
// needs to happen after init() of computes
// b/c a compute::init() can delete a fix, e.g. compute chunk/atom
list_init(INITIAL_INTEGRATE,n_initial_integrate,list_initial_integrate);
list_init(POST_INTEGRATE,n_post_integrate,list_post_integrate);
list_init(PRE_EXCHANGE,n_pre_exchange,list_pre_exchange);
list_init(PRE_NEIGHBOR,n_pre_neighbor,list_pre_neighbor);
list_init(POST_NEIGHBOR,n_post_neighbor,list_post_neighbor);
list_init(PRE_FORCE,n_pre_force,list_pre_force);
list_init(PRE_REVERSE,n_pre_reverse,list_pre_reverse);
list_init(POST_FORCE,n_post_force,list_post_force);
list_init(FINAL_INTEGRATE,n_final_integrate,list_final_integrate);
list_init_end_of_step(END_OF_STEP,n_end_of_step,list_end_of_step);
list_init_thermo_energy(THERMO_ENERGY,n_thermo_energy,list_thermo_energy);
list_init_thermo_energy_atom(n_thermo_energy_atom,list_thermo_energy_atom);
list_init(INITIAL_INTEGRATE_RESPA,
n_initial_integrate_respa,list_initial_integrate_respa);
list_init(POST_INTEGRATE_RESPA,
n_post_integrate_respa,list_post_integrate_respa);
list_init(POST_FORCE_RESPA,
n_post_force_respa,list_post_force_respa);
list_init(PRE_FORCE_RESPA,
n_pre_force_respa,list_pre_force_respa);
list_init(FINAL_INTEGRATE_RESPA,
n_final_integrate_respa,list_final_integrate_respa);
list_init(MIN_PRE_EXCHANGE,n_min_pre_exchange,list_min_pre_exchange);
list_init(MIN_PRE_NEIGHBOR,n_min_pre_neighbor,list_min_pre_neighbor);
list_init(MIN_POST_NEIGHBOR,n_min_post_neighbor,list_min_post_neighbor);
list_init(MIN_PRE_FORCE,n_min_pre_force,list_min_pre_force);
list_init(MIN_PRE_REVERSE,n_min_pre_reverse,list_min_pre_reverse);
list_init(MIN_POST_FORCE,n_min_post_force,list_min_post_force);
list_init(MIN_ENERGY,n_min_energy,list_min_energy);
// create list of computes that store invocation times
list_init_compute();
// error if any fix or compute is using a dynamic group when not allowed
for (i = 0; i < nfix; i++)
if (!fix[i]->dynamic_group_allow && group->dynamic[fix[i]->igroup]) {
char str[128];
snprintf(str,128,
"Fix %s does not allow use of dynamic group",fix[i]->id);
error->all(FLERR,str);
}
for (i = 0; i < ncompute; i++)
if (!compute[i]->dynamic_group_allow &&
group->dynamic[compute[i]->igroup]) {
char str[128];
snprintf(str,128,"Compute %s does not allow use of dynamic group",
fix[i]->id);
error->all(FLERR,str);
}
// warn if any particle is time integrated more than once
int nlocal = atom->nlocal;
int *mask = atom->mask;
int *flag = new int[nlocal];
for (i = 0; i < nlocal; i++) flag[i] = 0;
int groupbit;
for (i = 0; i < nfix; i++) {
if (fix[i]->time_integrate == 0) continue;
groupbit = fix[i]->groupbit;
for (j = 0; j < nlocal; j++)
if (mask[j] & groupbit) flag[j]++;
}
int check = 0;
for (i = 0; i < nlocal; i++)
if (flag[i] > 1) check = 1;
delete [] flag;
int checkall;
MPI_Allreduce(&check,&checkall,1,MPI_INT,MPI_SUM,world);
if (comm->me == 0 && checkall)
error->warning(FLERR,
"One or more atoms are time integrated more than once");
}
/* ----------------------------------------------------------------------
setup for run, calls setup() of all fixes and computes
called from Verlet, RESPA, Min
------------------------------------------------------------------------- */
void Modify::setup(int vflag)
{
// compute setup needs to come before fix setup
// b/c NH fixes need DOF of temperature computes
// fix group setup() is special case since populates a dynamic group
// needs to be done before temperature compute setup
for (int i = 0; i < nfix; i++)
if (strcmp(fix[i]->style,"GROUP") == 0) fix[i]->setup(vflag);
for (int i = 0; i < ncompute; i++) compute[i]->setup();
if (update->whichflag == 1)
for (int i = 0; i < nfix; i++) fix[i]->setup(vflag);
else if (update->whichflag == 2)
for (int i = 0; i < nfix; i++) fix[i]->min_setup(vflag);
}
/* ----------------------------------------------------------------------
setup pre_exchange call, only for fixes that define pre_exchange
called from Verlet, RESPA, Min, and WriteRestart with whichflag = 0
------------------------------------------------------------------------- */
void Modify::setup_pre_exchange()
{
if (update->whichflag <= 1)
for (int i = 0; i < n_pre_exchange; i++)
fix[list_pre_exchange[i]]->setup_pre_exchange();
else if (update->whichflag == 2)
for (int i = 0; i < n_min_pre_exchange; i++)
fix[list_min_pre_exchange[i]]->setup_pre_exchange();
}
/* ----------------------------------------------------------------------
setup pre_neighbor call, only for fixes that define pre_neighbor
called from Verlet, RESPA
------------------------------------------------------------------------- */
void Modify::setup_pre_neighbor()
{
if (update->whichflag == 1)
for (int i = 0; i < n_pre_neighbor; i++)
fix[list_pre_neighbor[i]]->setup_pre_neighbor();
else if (update->whichflag == 2)
for (int i = 0; i < n_min_pre_neighbor; i++)
fix[list_min_pre_neighbor[i]]->setup_pre_neighbor();
}
/* ----------------------------------------------------------------------
setup post_neighbor call, only for fixes that define post_neighbor
called from Verlet, RESPA
------------------------------------------------------------------------- */
void Modify::setup_post_neighbor()
{
if (update->whichflag == 1)
for (int i = 0; i < n_post_neighbor; i++)
fix[list_post_neighbor[i]]->setup_post_neighbor();
else if (update->whichflag == 2)
for (int i = 0; i < n_min_post_neighbor; i++)
fix[list_min_post_neighbor[i]]->setup_post_neighbor();
}
/* ----------------------------------------------------------------------
setup pre_force call, only for fixes that define pre_force
called from Verlet, RESPA, Min
------------------------------------------------------------------------- */
void Modify::setup_pre_force(int vflag)
{
if (update->whichflag == 1)
for (int i = 0; i < n_pre_force; i++)
fix[list_pre_force[i]]->setup_pre_force(vflag);
else if (update->whichflag == 2)
for (int i = 0; i < n_min_pre_force; i++)
fix[list_min_pre_force[i]]->setup_pre_force(vflag);
}
/* ----------------------------------------------------------------------
setup pre_reverse call, only for fixes that define pre_reverse
called from Verlet, RESPA, Min
------------------------------------------------------------------------- */
void Modify::setup_pre_reverse(int eflag, int vflag)
{
if (update->whichflag == 1)
for (int i = 0; i < n_pre_reverse; i++)
fix[list_pre_reverse[i]]->setup_pre_reverse(eflag,vflag);
else if (update->whichflag == 2)
for (int i = 0; i < n_min_pre_reverse; i++)
fix[list_min_pre_reverse[i]]->setup_pre_reverse(eflag,vflag);
}
/* ----------------------------------------------------------------------
1st half of integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::initial_integrate(int vflag)
{
for (int i = 0; i < n_initial_integrate; i++)
fix[list_initial_integrate[i]]->initial_integrate(vflag);
}
/* ----------------------------------------------------------------------
post_integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::post_integrate()
{
for (int i = 0; i < n_post_integrate; i++)
fix[list_post_integrate[i]]->post_integrate();
}
/* ----------------------------------------------------------------------
pre_exchange call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::pre_exchange()
{
for (int i = 0; i < n_pre_exchange; i++)
fix[list_pre_exchange[i]]->pre_exchange();
}
/* ----------------------------------------------------------------------
pre_neighbor call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::pre_neighbor()
{
for (int i = 0; i < n_pre_neighbor; i++)
fix[list_pre_neighbor[i]]->pre_neighbor();
}
/* ----------------------------------------------------------------------
post_neighbor call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::post_neighbor()
{
for (int i = 0; i < n_post_neighbor; i++)
fix[list_post_neighbor[i]]->post_neighbor();
}
/* ----------------------------------------------------------------------
pre_force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::pre_force(int vflag)
{
for (int i = 0; i < n_pre_force; i++)
fix[list_pre_force[i]]->pre_force(vflag);
}
/* ----------------------------------------------------------------------
pre_reverse call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::pre_reverse(int eflag, int vflag)
{
for (int i = 0; i < n_pre_reverse; i++)
fix[list_pre_reverse[i]]->pre_reverse(eflag,vflag);
}
/* ----------------------------------------------------------------------
post_force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::post_force(int vflag)
{
for (int i = 0; i < n_post_force; i++)
fix[list_post_force[i]]->post_force(vflag);
}
/* ----------------------------------------------------------------------
2nd half of integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::final_integrate()
{
for (int i = 0; i < n_final_integrate; i++)
fix[list_final_integrate[i]]->final_integrate();
}
/* ----------------------------------------------------------------------
end-of-timestep call, only for relevant fixes
only call fix->end_of_step() on timesteps that are multiples of nevery
------------------------------------------------------------------------- */
void Modify::end_of_step()
{
for (int i = 0; i < n_end_of_step; i++)
if (update->ntimestep % end_of_step_every[i] == 0)
fix[list_end_of_step[i]]->end_of_step();
}
/* ----------------------------------------------------------------------
thermo energy call, only for relevant fixes
called by Thermo class
compute_scalar() is fix call to return energy
------------------------------------------------------------------------- */
double Modify::thermo_energy()
{
double energy = 0.0;
for (int i = 0; i < n_thermo_energy; i++)
energy += fix[list_thermo_energy[i]]->compute_scalar();
return energy;
}
/* ----------------------------------------------------------------------
per-atom thermo energy call, only for relevant fixes
called by compute pe/atom
------------------------------------------------------------------------- */
void Modify::thermo_energy_atom(int nlocal, double *energy)
{
int i,j;
double *eatom;
for (i = 0; i < n_thermo_energy_atom; i++) {
eatom = fix[list_thermo_energy_atom[i]]->eatom;
if (!eatom) continue;
for (j = 0; j < nlocal; j++) energy[j] += eatom[j];
}
}
/* ----------------------------------------------------------------------
post_run call
------------------------------------------------------------------------- */
void Modify::post_run()
{
for (int i = 0; i < nfix; i++) fix[i]->post_run();
}
/* ----------------------------------------------------------------------
create_attribute call
invoked when an atom is added to system during a run
necessary so that fixes and computes that store per-atom
state can initialize that state for the new atom N
computes can store per-atom state via a fix like fix STORE
compute has the create_attribute flag, not fix STORE
------------------------------------------------------------------------- */
void Modify::create_attribute(int n)
{
for (int i = 0; i < nfix; i++)
if (fix[i]->create_attribute) fix[i]->set_arrays(n);
for (int i = 0; i < ncompute; i++)
if (compute[i]->create_attribute) compute[i]->set_arrays(n);
input->variable->set_arrays(n);
}
/* ----------------------------------------------------------------------
setup rRESPA pre_force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::setup_pre_force_respa(int vflag, int ilevel)
{
for (int i = 0; i < n_pre_force_respa; i++)
fix[list_pre_force_respa[i]]->setup_pre_force_respa(vflag,ilevel);
}
/* ----------------------------------------------------------------------
1st half of rRESPA integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::initial_integrate_respa(int vflag, int ilevel, int iloop)
{
for (int i = 0; i < n_initial_integrate_respa; i++)
fix[list_initial_integrate_respa[i]]->
initial_integrate_respa(vflag,ilevel,iloop);
}
/* ----------------------------------------------------------------------
rRESPA post_integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::post_integrate_respa(int ilevel, int iloop)
{
for (int i = 0; i < n_post_integrate_respa; i++)
fix[list_post_integrate_respa[i]]->post_integrate_respa(ilevel,iloop);
}
/* ----------------------------------------------------------------------
rRESPA pre_force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::pre_force_respa(int vflag, int ilevel, int iloop)
{
for (int i = 0; i < n_pre_force_respa; i++)
fix[list_pre_force_respa[i]]->pre_force_respa(vflag,ilevel,iloop);
}
/* ----------------------------------------------------------------------
rRESPA post_force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::post_force_respa(int vflag, int ilevel, int iloop)
{
for (int i = 0; i < n_post_force_respa; i++)
fix[list_post_force_respa[i]]->post_force_respa(vflag,ilevel,iloop);
}
/* ----------------------------------------------------------------------
2nd half of rRESPA integrate call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::final_integrate_respa(int ilevel, int iloop)
{
for (int i = 0; i < n_final_integrate_respa; i++)
fix[list_final_integrate_respa[i]]->final_integrate_respa(ilevel,iloop);
}
/* ----------------------------------------------------------------------
minimizer pre-exchange call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_pre_exchange()
{
for (int i = 0; i < n_min_pre_exchange; i++)
fix[list_min_pre_exchange[i]]->min_pre_exchange();
}
/* ----------------------------------------------------------------------
minimizer pre-neighbor call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_pre_neighbor()
{
for (int i = 0; i < n_min_pre_neighbor; i++)
fix[list_min_pre_neighbor[i]]->min_pre_neighbor();
}
/* ----------------------------------------------------------------------
minimizer post-neighbor call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_post_neighbor()
{
for (int i = 0; i < n_min_post_neighbor; i++)
fix[list_min_post_neighbor[i]]->min_post_neighbor();
}
/* ----------------------------------------------------------------------
minimizer pre-force call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_pre_force(int vflag)
{
for (int i = 0; i < n_min_pre_force; i++)
fix[list_min_pre_force[i]]->min_pre_force(vflag);
}
/* ----------------------------------------------------------------------
minimizer pre-reverse call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_pre_reverse(int eflag, int vflag)
{
for (int i = 0; i < n_min_pre_reverse; i++)
fix[list_min_pre_reverse[i]]->min_pre_reverse(eflag,vflag);
}
/* ----------------------------------------------------------------------
minimizer force adjustment call, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_post_force(int vflag)
{
for (int i = 0; i < n_min_post_force; i++)
fix[list_min_post_force[i]]->min_post_force(vflag);
}
/* ----------------------------------------------------------------------
minimizer energy/force evaluation, only for relevant fixes
return energy and forces on extra degrees of freedom
------------------------------------------------------------------------- */
double Modify::min_energy(double *fextra)
{
int ifix,index;
index = 0;
double eng = 0.0;
for (int i = 0; i < n_min_energy; i++) {
ifix = list_min_energy[i];
eng += fix[ifix]->min_energy(&fextra[index]);
index += fix[ifix]->min_dof();
}
return eng;
}
/* ----------------------------------------------------------------------
store current state of extra minimizer dof, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_store()
{
for (int i = 0; i < n_min_energy; i++)
fix[list_min_energy[i]]->min_store();
}
/* ----------------------------------------------------------------------
manage state of extra minimizer dof on a stack, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_clearstore()
{
for (int i = 0; i < n_min_energy; i++)
fix[list_min_energy[i]]->min_clearstore();
}
void Modify::min_pushstore()
{
for (int i = 0; i < n_min_energy; i++)
fix[list_min_energy[i]]->min_pushstore();
}
void Modify::min_popstore()
{
for (int i = 0; i < n_min_energy; i++)
fix[list_min_energy[i]]->min_popstore();
}
/* ----------------------------------------------------------------------
displace extra minimizer dof along vector hextra, only for relevant fixes
------------------------------------------------------------------------- */
void Modify::min_step(double alpha, double *hextra)
{
int ifix,index;
index = 0;
for (int i = 0; i < n_min_energy; i++) {
ifix = list_min_energy[i];
fix[ifix]->min_step(alpha,&hextra[index]);
index += fix[ifix]->min_dof();
}
}
/* ----------------------------------------------------------------------
compute max allowed step size along vector hextra, only for relevant fixes
------------------------------------------------------------------------- */
double Modify::max_alpha(double *hextra)
{
int ifix,index;
double alpha = BIG;
index = 0;
for (int i = 0; i < n_min_energy; i++) {
ifix = list_min_energy[i];
double alpha_one = fix[ifix]->max_alpha(&hextra[index]);
alpha = MIN(alpha,alpha_one);
index += fix[ifix]->min_dof();
}
return alpha;
}
/* ----------------------------------------------------------------------
extract extra minimizer dof, only for relevant fixes
------------------------------------------------------------------------- */
int Modify::min_dof()
{
int ndof = 0;
for (int i = 0; i < n_min_energy; i++)
ndof += fix[list_min_energy[i]]->min_dof();
return ndof;
}
/* ----------------------------------------------------------------------
reset minimizer reference state of fix, only for relevant fixes
------------------------------------------------------------------------- */
int Modify::min_reset_ref()
{
int itmp,itmpall;
itmpall = 0;
for (int i = 0; i < n_min_energy; i++) {
itmp = fix[list_min_energy[i]]->min_reset_ref();
if (itmp) itmpall = 1;
}
return itmpall;
}
/* ----------------------------------------------------------------------
add a new fix or replace one with same ID
------------------------------------------------------------------------- */
void Modify::add_fix(int narg, char **arg, int trysuffix)
{
if (narg < 3) error->all(FLERR,"Illegal fix command");
// cannot define fix before box exists unless style is in exception list
// don't like this way of checking for exceptions by adding fixes to list,
// but can't think of better way
// too late if instantiate fix, then check flag set in fix constructor,
// since some fixes access domain settings in their constructor
// NULL must be last entry in this list
const char *exceptions[] =
{"GPU", "OMP", "INTEL", "property/atom", "cmap", "cmap3", "rx",
"deprecated", "STORE/KIM", NULL};
if (domain->box_exist == 0) {
int m;
for (m = 0; exceptions[m] != NULL; m++)
if (strcmp(arg[2],exceptions[m]) == 0) break;
if (exceptions[m] == NULL)
error->all(FLERR,"Fix command before simulation box is defined");
}
// check group ID
int igroup = group->find(arg[1]);
if (igroup == -1) error->all(FLERR,"Could not find fix group ID");
// if fix ID exists:
// set newflag = 0 so create new fix in same location in fix list
// error if new style does not match old style
// since can't replace it (all when-to-invoke ptrs would be invalid)
// warn if new group != old group
// delete old fix, but do not call update_callback(),
// since will replace this fix and thus other fix locs will not change
// set ptr to NULL in case new fix scans list of fixes,
// e.g. scan will occur in add_callback() if called by new fix
// if fix ID does not exist:
// set newflag = 1 so create new fix
// extend fix and fmask lists as necessary
int ifix,newflag;
for (ifix = 0; ifix < nfix; ifix++)
if (strcmp(arg[0],fix[ifix]->id) == 0) break;
if (ifix < nfix) {
newflag = 0;
int match = 0;
if (strcmp(arg[2],fix[ifix]->style) == 0) match = 1;
if (!match && trysuffix && lmp->suffix_enable) {
char estyle[256];
if (lmp->suffix) {
sprintf(estyle,"%s/%s",arg[2],lmp->suffix);
if (strcmp(estyle,fix[ifix]->style) == 0) match = 1;
}
if (lmp->suffix2) {
sprintf(estyle,"%s/%s",arg[2],lmp->suffix2);
if (strcmp(estyle,fix[ifix]->style) == 0) match = 1;
}
}
if (!match) error->all(FLERR,
"Replacing a fix, but new style != old style");
if (fix[ifix]->igroup != igroup && comm->me == 0)
error->warning(FLERR,"Replacing a fix, but new group != old group");
delete fix[ifix];
fix[ifix] = NULL;
} else {
newflag = 1;
if (nfix == maxfix) {
maxfix += DELTA;
fix = (Fix **) memory->srealloc(fix,maxfix*sizeof(Fix *),"modify:fix");
memory->grow(fmask,maxfix,"modify:fmask");
}
}
// create the Fix
// try first with suffix appended
fix[ifix] = NULL;
if (trysuffix && lmp->suffix_enable) {
if (lmp->suffix) {
int n = strlen(arg[2])+strlen(lmp->suffix)+2;
char *estyle = new char[n];
sprintf(estyle,"%s/%s",arg[2],lmp->suffix);
if (fix_map->find(estyle) != fix_map->end()) {
FixCreator fix_creator = (*fix_map)[estyle];
fix[ifix] = fix_creator(lmp,narg,arg);
delete[] fix[ifix]->style;
fix[ifix]->style = estyle;
} else delete[] estyle;
}
if (fix[ifix] == NULL && lmp->suffix2) {
int n = strlen(arg[2])+strlen(lmp->suffix2)+2;
char *estyle = new char[n];
sprintf(estyle,"%s/%s",arg[2],lmp->suffix2);
if (fix_map->find(estyle) != fix_map->end()) {
FixCreator fix_creator = (*fix_map)[estyle];
fix[ifix] = fix_creator(lmp,narg,arg);
delete[] fix[ifix]->style;
fix[ifix]->style = estyle;
} else delete[] estyle;
}
}
if (fix[ifix] == NULL && fix_map->find(arg[2]) != fix_map->end()) {
FixCreator fix_creator = (*fix_map)[arg[2]];
fix[ifix] = fix_creator(lmp,narg,arg);
}
if (fix[ifix] == NULL)
error->all(FLERR,utils::check_packages_for_style("fix",arg[2],lmp).c_str());
// check if Fix is in restart_global list
// if yes, pass state info to the Fix so it can reset itself
for (int i = 0; i < nfix_restart_global; i++)
if (strcmp(id_restart_global[i],fix[ifix]->id) == 0 &&
strcmp(style_restart_global[i],fix[ifix]->style) == 0) {
fix[ifix]->restart(state_restart_global[i]);
used_restart_global[i] = 1;
if (comm->me == 0) {
if (screen)
fprintf(screen,"Resetting global fix info from restart file:\n");
if (logfile)
fprintf(logfile,"Resetting global fix info from restart file:\n");
if (screen) fprintf(screen," fix style: %s, fix ID: %s\n",
fix[ifix]->style,fix[ifix]->id);
if (logfile) fprintf(logfile," fix style: %s, fix ID: %s\n",
fix[ifix]->style,fix[ifix]->id);
}
}
// check if Fix is in restart_peratom list
// if yes, loop over atoms so they can extract info from atom->extra array
for (int i = 0; i < nfix_restart_peratom; i++)
if (strcmp(id_restart_peratom[i],fix[ifix]->id) == 0 &&
strcmp(style_restart_peratom[i],fix[ifix]->style) == 0) {
used_restart_peratom[i] = 1;
for (int j = 0; j < atom->nlocal; j++)
fix[ifix]->unpack_restart(j,index_restart_peratom[i]);
fix[ifix]->restart_reset = 1;
if (comm->me == 0) {
if (screen)
fprintf(screen,"Resetting peratom fix info from restart file:\n");
if (logfile)
fprintf(logfile,"Resetting peratom fix info from restart file:\n");
if (screen) fprintf(screen," fix style: %s, fix ID: %s\n",
fix[ifix]->style,fix[ifix]->id);
if (logfile) fprintf(logfile," fix style: %s, fix ID: %s\n",
fix[ifix]->style,fix[ifix]->id);
}
}
// increment nfix (if new)
// set fix mask values
// post_constructor() allows new fix to create other fixes
// nfix increment comes first so that recursive call to add_fix within
// post_constructor() will see updated nfix
if (newflag) nfix++;
fmask[ifix] = fix[ifix]->setmask();
fix[ifix]->post_constructor();
}
/* ----------------------------------------------------------------------
one instance per fix in style_fix.h
------------------------------------------------------------------------- */
template <typename T>
Fix *Modify::fix_creator(LAMMPS *lmp, int narg, char **arg)
{
return new T(lmp,narg,arg);
}
/* ----------------------------------------------------------------------
modify a Fix's parameters
------------------------------------------------------------------------- */
void Modify::modify_fix(int narg, char **arg)
{
if (narg < 2) error->all(FLERR,"Illegal fix_modify command");
// lookup Fix ID
int ifix;
for (ifix = 0; ifix < nfix; ifix++)
if (strcmp(arg[0],fix[ifix]->id) == 0) break;
if (ifix == nfix) error->all(FLERR,"Could not find fix_modify ID");
fix[ifix]->modify_params(narg-1,&arg[1]);
}
/* ----------------------------------------------------------------------
delete a Fix from list of Fixes
Atom class must update indices in its list of callbacks to fixes
------------------------------------------------------------------------- */
void Modify::delete_fix(const char *id)
{
int ifix = find_fix(id);
if (ifix < 0) error->all(FLERR,"Could not find fix ID to delete");
delete_fix(ifix);
}
void Modify::delete_fix(int ifix)
{
if (fix[ifix]) delete fix[ifix];
atom->update_callback(ifix);
// move other Fixes and fmask down in list one slot
for (int i = ifix+1; i < nfix; i++) fix[i-1] = fix[i];
for (int i = ifix+1; i < nfix; i++) fmask[i-1] = fmask[i];
nfix--;
}
/* ----------------------------------------------------------------------
find a fix by ID
return index of fix or -1 if not found
------------------------------------------------------------------------- */
int Modify::find_fix(const char *id)
{
if (id == NULL) return -1;
int ifix;
for (ifix = 0; ifix < nfix; ifix++)
if (strcmp(id,fix[ifix]->id) == 0) break;
if (ifix == nfix) return -1;
return ifix;
}
/* ----------------------------------------------------------------------
find a fix by style
return index of fix or -1 if not found
------------------------------------------------------------------------- */
int Modify::find_fix_by_style(const char *style)
{
int ifix;
for (ifix = 0; ifix < nfix; ifix++)
if (strcmp(style,fix[ifix]->style) == 0) break;
if (ifix == nfix) return -1;
return ifix;
}
/* ----------------------------------------------------------------------
check for fix associated with package name in compiled list
return 1 if found else 0
used to determine whether LAMMPS was built with
GPU, USER-INTEL, USER-OMP packages, which have their own fixes
------------------------------------------------------------------------- */
int Modify::check_package(const char *package_fix_name)
{
if (fix_map->find(package_fix_name) == fix_map->end()) return 0;
return 1;
}
/* ----------------------------------------------------------------------
check if the group indicated by groupbit overlaps with any
currently existing rigid fixes. return 1 in this case otherwise 0
------------------------------------------------------------------------- */
int Modify::check_rigid_group_overlap(int groupbit)
{
const int * const mask = atom->mask;
const int nlocal = atom->nlocal;
int dim;
int n = 0;
for (int ifix = 0; ifix < nfix; ifix++) {
if (utils::strmatch(fix[ifix]->style,"^rigid")) {
const int * const body = (const int *)fix[ifix]->extract("body",dim);
if ((body == NULL) || (dim != 1)) break;
for (int i=0; (i < nlocal) && (n == 0); ++i)
if ((mask[i] & groupbit) && (body[i] >= 0)) ++n;
}
}
int n_all = 0;
MPI_Allreduce(&n,&n_all,1,MPI_INT,MPI_SUM,world);
if (n_all > 0) return 1;
return 0;
}
/* ----------------------------------------------------------------------
check if the atoms in the group indicated by groupbit _and_ region
indicated by regionid overlap with any currently existing rigid fixes.
return 1 in this case, otherwise 0
------------------------------------------------------------------------- */
int Modify::check_rigid_region_overlap(int groupbit, Region *reg)
{
const int * const mask = atom->mask;
const double * const * const x = atom->x;
const int nlocal = atom->nlocal;
int dim;
int n = 0;
reg->prematch();
for (int ifix = 0; ifix < nfix; ifix++) {
if (strncmp("rigid",fix[ifix]->style,5) == 0) {
const int * const body = (const int *)fix[ifix]->extract("body",dim);
if ((body == NULL) || (dim != 1)) break;
for (int i=0; (i < nlocal) && (n == 0); ++i)
if ((mask[i] & groupbit) && (body[i] >= 0)
&& reg->match(x[i][0],x[i][1],x[i][2])) ++n;
}
}
int n_all = 0;
MPI_Allreduce(&n,&n_all,1,MPI_INT,MPI_SUM,world);
if (n_all > 0) return 1;
return 0;
}
/* ----------------------------------------------------------------------
check if the atoms in the selection list (length atom->nlocal,
content: 1 if atom is contained, 0 if not) overlap with currently
existing rigid fixes. return 1 in this case otherwise 0
------------------------------------------------------------------------- */
int Modify::check_rigid_list_overlap(int *select)
{
const int nlocal = atom->nlocal;
int dim;
int n = 0;
for (int ifix = 0; ifix < nfix; ifix++) {
if (strncmp("rigid",fix[ifix]->style,5) == 0) {
const int * const body = (const int *)fix[ifix]->extract("body",dim);
if ((body == NULL) || (dim != 1)) break;
for (int i=0; (i < nlocal) && (n == 0); ++i)
if ((body[i] >= 0) && select[i]) ++n;
}
}
int n_all = 0;
MPI_Allreduce(&n,&n_all,1,MPI_INT,MPI_SUM,world);
if (n_all > 0) return 1;
return 0;
}
/* ----------------------------------------------------------------------
add a new compute
------------------------------------------------------------------------- */
void Modify::add_compute(int narg, char **arg, int trysuffix)
{
if (narg < 3) error->all(FLERR,"Illegal compute command");
// error check
for (int icompute = 0; icompute < ncompute; icompute++)
if (strcmp(arg[0],compute[icompute]->id) == 0)
error->all(FLERR,"Reuse of compute ID");
// extend Compute list if necessary
if (ncompute == maxcompute) {
maxcompute += DELTA;
compute = (Compute **)
memory->srealloc(compute,maxcompute*sizeof(Compute *),"modify:compute");
}
// create the Compute
// try first with suffix appended
compute[ncompute] = NULL;
if (trysuffix && lmp->suffix_enable) {
if (lmp->suffix) {
int n = strlen(arg[2])+strlen(lmp->suffix)+2;
char *estyle = new char[n];
sprintf(estyle,"%s/%s",arg[2],lmp->suffix);
if (compute_map->find(estyle) != compute_map->end()) {
ComputeCreator compute_creator = (*compute_map)[estyle];
compute[ncompute] = compute_creator(lmp,narg,arg);
delete[] compute[ncompute]->style;
compute[ncompute]->style = estyle;
} else delete[] estyle;
}
if (compute[ncompute] == NULL && lmp->suffix2) {
int n = strlen(arg[2])+strlen(lmp->suffix2)+2;
char *estyle = new char[n];
sprintf(estyle,"%s/%s",arg[2],lmp->suffix2);
if (compute_map->find(estyle) != compute_map->end()) {
ComputeCreator compute_creator = (*compute_map)[estyle];
compute[ncompute] = compute_creator(lmp,narg,arg);
delete[] compute[ncompute]->style;
compute[ncompute]->style = estyle;
} else delete[] estyle;
}
}
if (compute[ncompute] == NULL &&
compute_map->find(arg[2]) != compute_map->end()) {
ComputeCreator compute_creator = (*compute_map)[arg[2]];
compute[ncompute] = compute_creator(lmp,narg,arg);
}
if (compute[ncompute] == NULL)
error->all(FLERR,utils::check_packages_for_style("compute",arg[2],lmp).c_str());
ncompute++;
}
/* ----------------------------------------------------------------------
one instance per compute in style_compute.h
------------------------------------------------------------------------- */
template <typename T>
Compute *Modify::compute_creator(LAMMPS *lmp, int narg, char **arg)
{
return new T(lmp,narg,arg);
}
/* ----------------------------------------------------------------------
modify a Compute's parameters
------------------------------------------------------------------------- */
void Modify::modify_compute(int narg, char **arg)
{
if (narg < 2) error->all(FLERR,"Illegal compute_modify command");
// lookup Compute ID
int icompute;
for (icompute = 0; icompute < ncompute; icompute++)
if (strcmp(arg[0],compute[icompute]->id) == 0) break;
if (icompute == ncompute)
error->all(FLERR,"Could not find compute_modify ID");
compute[icompute]->modify_params(narg-1,&arg[1]);
}
/* ----------------------------------------------------------------------
delete a Compute from list of Computes
------------------------------------------------------------------------- */
void Modify::delete_compute(const char *id)
{
int icompute = find_compute(id);
if (icompute < 0) error->all(FLERR,"Could not find compute ID to delete");
delete compute[icompute];
// move other Computes down in list one slot
for (int i = icompute+1; i < ncompute; i++) compute[i-1] = compute[i];
ncompute--;
}
/* ----------------------------------------------------------------------
find a compute by ID
return index of compute or -1 if not found
------------------------------------------------------------------------- */
int Modify::find_compute(const char *id)
{
if(id==NULL) return -1;
int icompute;
for (icompute = 0; icompute < ncompute; icompute++)
if (strcmp(id,compute[icompute]->id) == 0) break;
if (icompute == ncompute) return -1;
return icompute;
}
/* ----------------------------------------------------------------------
clear invoked flag of all computes
called everywhere that computes are used, before computes are invoked
invoked flag used to avoid re-invoking same compute multiple times
and to flag computes that store invocation times as having been invoked
------------------------------------------------------------------------- */
void Modify::clearstep_compute()
{
for (int icompute = 0; icompute < ncompute; icompute++)
compute[icompute]->invoked_flag = 0;
}
/* ----------------------------------------------------------------------
loop over computes that store invocation times
if its invoked flag set on this timestep, schedule next invocation
called everywhere that computes are used, after computes are invoked
------------------------------------------------------------------------- */
void Modify::addstep_compute(bigint newstep)
{
for (int icompute = 0; icompute < n_timeflag; icompute++)
if (compute[list_timeflag[icompute]]->invoked_flag)
compute[list_timeflag[icompute]]->addstep(newstep);
}
/* ----------------------------------------------------------------------
loop over all computes
schedule next invocation for those that store invocation times
called when not sure what computes will be needed on newstep
do not loop only over n_timeflag, since may not be set yet
------------------------------------------------------------------------- */
void Modify::addstep_compute_all(bigint newstep)
{
for (int icompute = 0; icompute < ncompute; icompute++)
if (compute[icompute]->timeflag) compute[icompute]->addstep(newstep);
}
/* ----------------------------------------------------------------------
write to restart file for all Fixes with restart info
(1) fixes that have global state
(2) fixes that store per-atom quantities
------------------------------------------------------------------------- */
void Modify::write_restart(FILE *fp)
{
int me = comm->me;
int count = 0;
for (int i = 0; i < nfix; i++)
if (fix[i]->restart_global) count++;
if (me == 0) fwrite(&count,sizeof(int),1,fp);
int n;
for (int i = 0; i < nfix; i++)
if (fix[i]->restart_global) {
if (me == 0) {
n = strlen(fix[i]->id) + 1;
fwrite(&n,sizeof(int),1,fp);
fwrite(fix[i]->id,sizeof(char),n,fp);
n = strlen(fix[i]->style) + 1;
fwrite(&n,sizeof(int),1,fp);
fwrite(fix[i]->style,sizeof(char),n,fp);
}
fix[i]->write_restart(fp);
}
count = 0;
for (int i = 0; i < nfix; i++)
if (fix[i]->restart_peratom) count++;
if (me == 0) fwrite(&count,sizeof(int),1,fp);
for (int i = 0; i < nfix; i++)
if (fix[i]->restart_peratom) {
int maxsize_restart = fix[i]->maxsize_restart();
if (me == 0) {
n = strlen(fix[i]->id) + 1;
fwrite(&n,sizeof(int),1,fp);
fwrite(fix[i]->id,sizeof(char),n,fp);
n = strlen(fix[i]->style) + 1;
fwrite(&n,sizeof(int),1,fp);
fwrite(fix[i]->style,sizeof(char),n,fp);
fwrite(&maxsize_restart,sizeof(int),1,fp);
}
}
}
/* ----------------------------------------------------------------------
read in restart file data on all previously defined Fixes with restart info
(1) fixes that have global state
(2) fixes that store per-atom quantities
return maxsize of extra info that will be stored with any atom
------------------------------------------------------------------------- */
int Modify::read_restart(FILE *fp)
{
// nfix_restart_global = # of restart entries with global state info
int me = comm->me;
if (me == 0) fread(&nfix_restart_global,sizeof(int),1,fp);
MPI_Bcast(&nfix_restart_global,1,MPI_INT,0,world);
// allocate space for each entry
if (nfix_restart_global) {
id_restart_global = new char*[nfix_restart_global];
style_restart_global = new char*[nfix_restart_global];
state_restart_global = new char*[nfix_restart_global];
used_restart_global = new int[nfix_restart_global];
}
// read each entry and Bcast to all procs
// each entry has id string, style string, chunk of state data
int n;
for (int i = 0; i < nfix_restart_global; i++) {
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
id_restart_global[i] = new char[n];
if (me == 0) fread(id_restart_global[i],sizeof(char),n,fp);
MPI_Bcast(id_restart_global[i],n,MPI_CHAR,0,world);
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
style_restart_global[i] = new char[n];
if (me == 0) fread(style_restart_global[i],sizeof(char),n,fp);
MPI_Bcast(style_restart_global[i],n,MPI_CHAR,0,world);
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
state_restart_global[i] = new char[n];
if (me == 0) fread(state_restart_global[i],sizeof(char),n,fp);
MPI_Bcast(state_restart_global[i],n,MPI_CHAR,0,world);
used_restart_global[i] = 0;
}
// nfix_restart_peratom = # of restart entries with peratom info
int maxsize = 0;
if (me == 0) fread(&nfix_restart_peratom,sizeof(int),1,fp);
MPI_Bcast(&nfix_restart_peratom,1,MPI_INT,0,world);
// allocate space for each entry
if (nfix_restart_peratom) {
id_restart_peratom = new char*[nfix_restart_peratom];
style_restart_peratom = new char*[nfix_restart_peratom];
index_restart_peratom = new int[nfix_restart_peratom];
used_restart_peratom = new int[nfix_restart_peratom];
}
// read each entry and Bcast to all procs
// each entry has id string, style string, maxsize of one atom's data
// set index = which set of extra data this fix represents
for (int i = 0; i < nfix_restart_peratom; i++) {
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
id_restart_peratom[i] = new char[n];
if (me == 0) fread(id_restart_peratom[i],sizeof(char),n,fp);
MPI_Bcast(id_restart_peratom[i],n,MPI_CHAR,0,world);
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
style_restart_peratom[i] = new char[n];
if (me == 0) fread(style_restart_peratom[i],sizeof(char),n,fp);
MPI_Bcast(style_restart_peratom[i],n,MPI_CHAR,0,world);
if (me == 0) fread(&n,sizeof(int),1,fp);
MPI_Bcast(&n,1,MPI_INT,0,world);
maxsize += n;
index_restart_peratom[i] = i;
used_restart_peratom[i] = 0;
}
return maxsize;
}
/* ----------------------------------------------------------------------
delete all lists of restart file Fix info
if flag set, print list of restart file info not assigned to new fixes
------------------------------------------------------------------------- */
void Modify::restart_deallocate(int flag)
{
if (nfix_restart_global) {
if (flag && comm->me == 0) {
int i;
for (i = 0; i < nfix_restart_global; i++)
if (used_restart_global[i] == 0) break;
if (i == nfix_restart_global) {
if (screen)
fprintf(screen,"All restart file global fix info "
"was re-assigned\n");
if (logfile)
fprintf(logfile,"All restart file global fix info "
"was re-assigned\n");
} else {
if (screen) fprintf(screen,"Unused restart file global fix info:\n");
if (logfile) fprintf(logfile,"Unused restart file global fix info:\n");
for (i = 0; i < nfix_restart_global; i++) {
if (used_restart_global[i]) continue;
if (screen) fprintf(screen," fix style: %s, fix ID: %s\n",
style_restart_global[i],id_restart_global[i]);
if (logfile) fprintf(logfile," fix style: %s, fix ID: %s\n",
style_restart_global[i],id_restart_global[i]);
}
}
}
for (int i = 0; i < nfix_restart_global; i++) {
delete [] id_restart_global[i];
delete [] style_restart_global[i];
delete [] state_restart_global[i];
}
delete [] id_restart_global;
delete [] style_restart_global;
delete [] state_restart_global;
delete [] used_restart_global;
}
if (nfix_restart_peratom) {
if (flag && comm->me == 0) {
int i;
for (i = 0; i < nfix_restart_peratom; i++)
if (used_restart_peratom[i] == 0) break;
if (i == nfix_restart_peratom) {
if (screen)
fprintf(screen,"All restart file peratom fix info "
"was re-assigned\n");
if (logfile)
fprintf(logfile,"All restart file peratom fix info "
"was re-assigned\n");
} else {
if (screen) fprintf(screen,"Unused restart file peratom fix info:\n");
if (logfile) fprintf(logfile,"Unused restart file peratom fix info:\n");
for (i = 0; i < nfix_restart_peratom; i++) {
if (used_restart_peratom[i]) continue;
if (screen) fprintf(screen," fix style: %s, fix ID: %s\n",
style_restart_peratom[i],id_restart_peratom[i]);
if (logfile) fprintf(logfile," fix style: %s, fix ID: %s\n",
style_restart_peratom[i],id_restart_peratom[i]);
}
}
}
for (int i = 0; i < nfix_restart_peratom; i++) {
delete [] id_restart_peratom[i];
delete [] style_restart_peratom[i];
}
delete [] id_restart_peratom;
delete [] style_restart_peratom;
delete [] index_restart_peratom;
delete [] used_restart_peratom;
}
nfix_restart_global = nfix_restart_peratom = 0;
}
/* ----------------------------------------------------------------------
create list of fix indices for fixes which match mask
------------------------------------------------------------------------- */
void Modify::list_init(int mask, int &n, int *&list)
{
delete [] list;
n = 0;
for (int i = 0; i < nfix; i++) if (fmask[i] & mask) n++;
list = new int[n];
n = 0;
for (int i = 0; i < nfix; i++) if (fmask[i] & mask) list[n++] = i;
}
/* ----------------------------------------------------------------------
create list of fix indices for end_of_step fixes
also create end_of_step_every[]
------------------------------------------------------------------------- */
void Modify::list_init_end_of_step(int mask, int &n, int *&list)
{
delete [] list;
delete [] end_of_step_every;
n = 0;
for (int i = 0; i < nfix; i++) if (fmask[i] & mask) n++;
list = new int[n];
end_of_step_every = new int[n];
n = 0;
for (int i = 0; i < nfix; i++)
if (fmask[i] & mask) {
list[n] = i;
end_of_step_every[n++] = fix[i]->nevery;
}
}
/* ----------------------------------------------------------------------
create list of fix indices for thermo energy fixes
only added to list if fix has THERMO_ENERGY mask set,
and its thermo_energy flag was set via fix_modify
------------------------------------------------------------------------- */
void Modify::list_init_thermo_energy(int mask, int &n, int *&list)
{
delete [] list;
n = 0;
for (int i = 0; i < nfix; i++)
if (fmask[i] & mask && fix[i]->thermo_energy) n++;
list = new int[n];
n = 0;
for (int i = 0; i < nfix; i++)
if (fmask[i] & mask && fix[i]->thermo_energy) list[n++] = i;
}
/* ----------------------------------------------------------------------
create list of fix indices for peratom thermo energy fixes
only added to list if fix has its peatom_flag set,
and its thermo_energy flag was set via fix_modify
------------------------------------------------------------------------- */
void Modify::list_init_thermo_energy_atom(int &n, int *&list)
{
delete [] list;
n = 0;
for (int i = 0; i < nfix; i++)
if (fix[i]->peatom_flag && fix[i]->thermo_energy) n++;
list = new int[n];
n = 0;
for (int i = 0; i < nfix; i++)
if (fix[i]->peatom_flag && fix[i]->thermo_energy) list[n++] = i;
}
/* ----------------------------------------------------------------------
create list of compute indices for computes which store invocation times
------------------------------------------------------------------------- */
void Modify::list_init_compute()
{
delete [] list_timeflag;
n_timeflag = 0;
for (int i = 0; i < ncompute; i++)
if (compute[i]->timeflag) n_timeflag++;
list_timeflag = new int[n_timeflag];
n_timeflag = 0;
for (int i = 0; i < ncompute; i++)
if (compute[i]->timeflag) list_timeflag[n_timeflag++] = i;
}
/* ----------------------------------------------------------------------
return # of bytes of allocated memory from all fixes
------------------------------------------------------------------------- */
bigint Modify::memory_usage()
{
bigint bytes = 0;
for (int i = 0; i < nfix; i++)
bytes += static_cast<bigint> (fix[i]->memory_usage());
for (int i = 0; i < ncompute; i++)
bytes += static_cast<bigint> (compute[i]->memory_usage());
return bytes;
}