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 package ffx.potential;
39
40 import ffx.crystal.Crystal;
41 import ffx.crystal.SpaceGroup;
42 import ffx.numerics.Potential;
43 import ffx.potential.MolecularAssembly.FractionalMode;
44 import ffx.potential.bonded.Atom;
45
46 import java.util.ArrayList;
47 import java.util.List;
48
49 import static org.apache.commons.math3.util.FastMath.toDegrees;
50
51
52
53
54
55
56
57
58 public class XtalEnergy implements Potential {
59
60
61
62
63 private final MolecularAssembly molecularAssembly;
64
65
66
67 private final ForceFieldEnergy forceFieldEnergy;
68
69
70
71 private final Atom[] activeAtoms;
72
73
74
75 private final int nActive;
76
77
78
79 private final double[] xyz;
80 private final double[] gr;
81 private final int nParams;
82
83 private final Crystal crystal;
84 private final boolean latticeOnly;
85 private final VARIABLE_TYPE[] type;
86 private final double[] mass;
87 private final Crystal unitCell;
88 private double[] scaling;
89 private double totalEnergy;
90 private FractionalMode fractionalMode = FractionalMode.OFF;
91
92
93
94
95
96
97
98 public XtalEnergy(ForceFieldEnergy forceFieldEnergy, MolecularAssembly molecularAssembly) {
99 this(forceFieldEnergy, molecularAssembly, false);
100 }
101
102
103
104
105
106
107
108
109 public XtalEnergy(ForceFieldEnergy forceFieldEnergy, MolecularAssembly molecularAssembly, boolean latticeOnly) {
110 this.forceFieldEnergy = forceFieldEnergy;
111 this.molecularAssembly = molecularAssembly;
112 this.latticeOnly = latticeOnly;
113
114 if (!latticeOnly) {
115 Atom[] atoms = molecularAssembly.getAtomArray();
116 List<Atom> active = new ArrayList<>();
117 for (Atom a : atoms) {
118 if (a.isActive()) {
119 active.add(a);
120 }
121 }
122 nActive = active.size();
123 activeAtoms = active.toArray(new Atom[0]);
124 xyz = new double[3 * nActive];
125 gr = new double[3 * nActive];
126 } else {
127 nActive = 0;
128 activeAtoms = null;
129 xyz = null;
130 gr = null;
131 }
132
133 nParams = 3 * nActive + 6;
134 crystal = forceFieldEnergy.getCrystal();
135 unitCell = crystal.getUnitCell();
136 type = new VARIABLE_TYPE[nParams];
137 mass = new double[nParams];
138
139
140 int index = 0;
141 for (int i = 0; i < nActive; i++) {
142 double m = activeAtoms[i].getMass();
143 mass[index] = m;
144 mass[index + 1] = m;
145 mass[index + 2] = m;
146 type[index] = VARIABLE_TYPE.X;
147 type[index + 1] = VARIABLE_TYPE.Y;
148 type[index + 2] = VARIABLE_TYPE.Z;
149 index += 3;
150 }
151
152 for (int i = nActive * 3; i < nActive * 3 + 6; i++) {
153 mass[i] = 1.0;
154 type[i] = VARIABLE_TYPE.OTHER;
155 }
156 }
157
158
159
160
161 @Override
162 public boolean destroy() {
163 return forceFieldEnergy.destroy();
164 }
165
166
167
168
169 @Override
170 public double energy(double[] x) {
171
172
173 unscaleCoordinates(x);
174
175
176 setCoordinates(x);
177
178
179 totalEnergy = forceFieldEnergy.energy(false, false);
180
181
182 scaleCoordinates(x);
183
184 return totalEnergy;
185 }
186
187
188
189
190 @Override
191 public double energyAndGradient(double[] x, double[] g) {
192
193 unscaleCoordinates(x);
194
195
196 setCoordinates(x);
197
198 if (latticeOnly) {
199
200 totalEnergy = forceFieldEnergy.energy(false, false);
201 } else {
202
203 totalEnergy = forceFieldEnergy.energyAndGradient(xyz, gr);
204
205 }
206
207
208 packGradient(x, g);
209
210
211 unitCellParameterDerivatives(x, g);
212
213 return totalEnergy;
214 }
215
216
217
218
219 @Override
220 public double[] getAcceleration(double[] acceleration) {
221 throw new UnsupportedOperationException("Not supported yet.");
222 }
223
224
225
226
227 @Override
228 public double[] getCoordinates(double[] x) {
229 int n = getNumberOfVariables();
230 if (x == null || x.length < n) {
231 x = new double[n];
232 }
233 int index = 0;
234 for (int i = 0; i < nActive; i++) {
235 Atom a = activeAtoms[i];
236 x[index] = a.getX();
237 index++;
238 x[index] = a.getY();
239 index++;
240 x[index] = a.getZ();
241 index++;
242 }
243 x[index] = unitCell.a;
244 index++;
245 x[index] = unitCell.b;
246 index++;
247 x[index] = unitCell.c;
248 index++;
249 x[index] = unitCell.alpha;
250 index++;
251 x[index] = unitCell.beta;
252 index++;
253 x[index] = unitCell.gamma;
254 return x;
255 }
256
257
258
259
260 @Override
261 public STATE getEnergyTermState() {
262 return forceFieldEnergy.getEnergyTermState();
263 }
264
265
266
267
268 @Override
269 public void setEnergyTermState(STATE state) {
270 forceFieldEnergy.setEnergyTermState(state);
271 }
272
273
274
275
276 @Override
277 public double[] getMass() {
278 return mass;
279 }
280
281
282
283
284 @Override
285 public int getNumberOfVariables() {
286 return nParams;
287 }
288
289
290
291
292 @Override
293 public double[] getPreviousAcceleration(double[] previousAcceleration) {
294 throw new UnsupportedOperationException("Not supported yet.");
295 }
296
297
298
299
300 @Override
301 public double[] getScaling() {
302 return scaling;
303 }
304
305
306
307
308 @Override
309 public void setScaling(double[] scaling) {
310 this.scaling = scaling;
311 }
312
313
314
315
316 @Override
317 public double getTotalEnergy() {
318 return totalEnergy;
319 }
320
321
322
323
324 @Override
325 public VARIABLE_TYPE[] getVariableTypes() {
326 return type;
327 }
328
329
330
331
332 @Override
333 public double[] getVelocity(double[] velocity) {
334 throw new UnsupportedOperationException("Not supported yet.");
335 }
336
337
338
339
340 @Override
341 public void setAcceleration(double[] acceleration) {
342 throw new UnsupportedOperationException("Not supported yet.");
343 }
344
345
346
347
348
349
350 public void setFractionalCoordinateMode(FractionalMode fractionalMode) {
351 this.fractionalMode = fractionalMode;
352 molecularAssembly.setFractionalMode(fractionalMode);
353 }
354
355
356
357
358 @Override
359 public void setPreviousAcceleration(double[] previousAcceleration) {
360 throw new UnsupportedOperationException("Not supported yet.");
361 }
362
363
364
365
366 @Override
367 public void setVelocity(double[] velocity) {
368 throw new UnsupportedOperationException("Not supported yet.");
369 }
370
371
372
373
374
375
376
377 private void unitCellParameterDerivatives(double[] x, double[] g) {
378 int index = 3 * nActive;
379 double eps = 1.0e-5;
380 double deps = toDegrees(eps);
381 SpaceGroup spaceGroup = crystal.spaceGroup;
382 switch (spaceGroup.latticeSystem) {
383 case TRICLINIC_LATTICE -> {
384 g[index] = finiteDifference(x, index, eps);
385 index++;
386 g[index] = finiteDifference(x, index, eps);
387 index++;
388 g[index] = finiteDifference(x, index, eps);
389 index++;
390 g[index] = finiteDifference(x, index, deps);
391 index++;
392 g[index] = finiteDifference(x, index, deps);
393 index++;
394 g[index] = finiteDifference(x, index, deps);
395 }
396 case MONOCLINIC_LATTICE -> {
397
398 g[index] = finiteDifference(x, index, eps);
399 index++;
400 g[index] = finiteDifference(x, index, eps);
401 index++;
402 g[index] = finiteDifference(x, index, eps);
403 index++;
404 g[index] = 0.0;
405 index++;
406 g[index] = finiteDifference(x, index, deps);
407 index++;
408 g[index] = 0.0;
409 }
410 case ORTHORHOMBIC_LATTICE -> {
411
412 g[index] = finiteDifference(x, index, eps);
413 index++;
414 g[index] = finiteDifference(x, index, eps);
415 index++;
416 g[index] = finiteDifference(x, index, eps);
417 index++;
418 g[index] = 0.0;
419 index++;
420 g[index] = 0.0;
421 index++;
422 g[index] = 0.0;
423 }
424 case TETRAGONAL_LATTICE,
425 HEXAGONAL_LATTICE
426 -> {
427 g[index] = finiteDifference2(x, index, index + 1, eps);
428 index++;
429 g[index] = g[index - 1];
430 index++;
431 g[index] = finiteDifference(x, index, eps);
432 index++;
433 g[index] = 0.0;
434 index++;
435 g[index] = 0.0;
436 index++;
437 g[index] = 0.0;
438 }
439 case RHOMBOHEDRAL_LATTICE -> {
440
441 g[index] = finiteDifference3(x, index, index + 1, index + 2, eps);
442 index++;
443 g[index] = g[index - 1];
444 index++;
445 g[index] = g[index - 2];
446 index++;
447 g[index] = finiteDifference3(x, index, index + 1, index + 2, deps);
448 index++;
449 g[index] = g[index - 1];
450 index++;
451 g[index] = g[index - 2];
452 }
453
454 case CUBIC_LATTICE -> {
455
456 g[index] = finiteDifference3(x, index, index + 1, index + 2, eps);
457 index++;
458 g[index] = g[index - 1];
459 index++;
460 g[index] = g[index - 2];
461 index++;
462 g[index] = 0.0;
463 index++;
464 g[index] = 0.0;
465 index++;
466 g[index] = 0.0;
467 }
468 }
469
470
471
472 if (scaling != null) {
473 index = 3 * nActive;
474 for (int i = 0; i < 6; i++) {
475 g[index] /= scaling[index];
476 index++;
477 }
478 }
479 }
480
481
482
483
484
485
486
487
488
489 private double finiteDifference(double[] x, int index, double eps) {
490 double scale = 1.0;
491 if (scaling != null) {
492 scale = scaling[index];
493 }
494 double x1 = x[index];
495 final double param = x1 / scale;
496
497 x[index] = (param + eps / 2.0) * scale;
498 double ePlus = energy(x);
499 x[index] = (param - eps / 2.0) * scale;
500 double eMinus = energy(x);
501
502 x[index] = x1;
503
504 return (ePlus - eMinus) / eps;
505 }
506
507
508
509
510
511
512
513
514 private double finiteDifference2(double[] x, int index1, int index2, double eps) {
515 double scale1 = 1.0;
516 double scale2 = 1.0;
517
518 if (scaling != null) {
519 scale1 = scaling[index1];
520 scale2 = scaling[index2];
521 }
522
523 final double x1 = x[index1];
524 final double x2 = x[index2];
525 final double param1 = x1 / scale1;
526 final double param2 = x2 / scale2;
527 x[index1] = (param1 + eps / 2.0) * scale1;
528 x[index2] = (param2 + eps / 2.0) * scale2;
529 double ePlus = energy(x);
530 x[index1] = (param1 - eps / 2.0) * scale1;
531 x[index2] = (param2 - eps / 2.0) * scale2;
532 double eMinus = energy(x);
533
534 x[index1] = x1;
535 x[index2] = x2;
536
537 return (ePlus - eMinus) / eps;
538 }
539
540
541
542
543
544
545
546
547
548 private double finiteDifference3(double[] x, int index1, int index2, int index3, double eps) {
549 double scale1 = 1.0;
550 double scale2 = 1.0;
551 double scale3 = 1.0;
552 if (scaling != null) {
553 scale1 = scaling[index1];
554 scale2 = scaling[index2];
555 scale3 = scaling[index3];
556 }
557
558 final double x1 = x[index1];
559 final double x2 = x[index2];
560 final double x3 = x[index3];
561
562 final double param1 = x1 / scale1;
563 final double param2 = x2 / scale2;
564 final double param3 = x3 / scale3;
565 x[index1] = (param1 + eps / 2.0) * scale1;
566 x[index2] = (param2 + eps / 2.0) * scale2;
567 x[index3] = (param3 + eps / 2.0) * scale3;
568 final double ePlus = energy(x);
569 x[index1] = (param1 - eps / 2.0) * scale1;
570 x[index2] = (param2 - eps / 2.0) * scale2;
571 x[index3] = (param3 - eps / 2.0) * scale3;
572 final double eMinus = energy(x);
573
574 x[index1] = x1;
575 x[index2] = x2;
576 x[index3] = x3;
577
578 return (ePlus - eMinus) / eps;
579 }
580
581
582
583
584
585
586
587 private void packGradient(double[] x, double[] g) {
588
589 if (scaling != null) {
590 int len = x.length;
591 for (int i = 0; i < len; i++) {
592 if (fractionalMode == FractionalMode.OFF) {
593 g[i] /= scaling[i];
594 } else {
595
596
597 g[i] = 0.0;
598 }
599 x[i] *= scaling[i];
600 }
601 }
602 }
603
604
605
606
607
608
609 @Override
610 public void setCoordinates(double[] x) {
611 assert (x != null);
612
613
614 if (fractionalMode != FractionalMode.OFF) {
615 molecularAssembly.computeFractionalCoordinates();
616 }
617
618 int index = nActive * 3;
619 double a = x[index];
620 double b = x[index + 1];
621 double c = x[index + 2];
622 double alpha = x[index + 3];
623 double beta = x[index + 4];
624 double gamma = x[index + 5];
625
626 SpaceGroup spaceGroup = crystal.spaceGroup;
627
628
629 switch (spaceGroup.latticeSystem) {
630 case TRICLINIC_LATTICE -> {
631 }
632 case MONOCLINIC_LATTICE -> {
633
634 alpha = 90.0;
635 gamma = 90.0;
636 }
637 case ORTHORHOMBIC_LATTICE -> {
638
639 alpha = 90.0;
640 beta = 90.0;
641 gamma = 90.0;
642 }
643 case TETRAGONAL_LATTICE -> {
644
645 double ab = 0.5 * (a + b);
646 a = ab;
647 b = ab;
648 alpha = 90.0;
649 beta = 90.0;
650 gamma = 90.0;
651 }
652 case RHOMBOHEDRAL_LATTICE -> {
653
654 double abc = (a + b + c) / 3.0;
655 a = abc;
656 b = abc;
657 c = abc;
658 double angles = (alpha + beta + gamma) / 3.0;
659 alpha = angles;
660 beta = angles;
661 gamma = angles;
662 }
663 case HEXAGONAL_LATTICE -> {
664
665 double ab = 0.5 * (a + b);
666 a = ab;
667 b = ab;
668 alpha = 90.0;
669 beta = 90.0;
670 gamma = 120.0;
671 }
672 case CUBIC_LATTICE -> {
673
674 double abc = (a + b + c) / 3.0;
675 a = abc;
676 b = abc;
677 c = abc;
678 alpha = 90.0;
679 beta = 90.0;
680 gamma = 90.0;
681 }
682 }
683
684 crystal.changeUnitCellParameters(a, b, c, alpha, beta, gamma);
685 forceFieldEnergy.setCrystal(crystal);
686
687
688 if (fractionalMode == FractionalMode.OFF) {
689 index = 0;
690 for (int i = 0; i < nActive; i++) {
691 Atom atom = activeAtoms[i];
692 double xx = x[index];
693 double yy = x[index + 1];
694 double zz = x[index + 2];
695 xyz[index] = xx;
696 xyz[index + 1] = yy;
697 xyz[index + 2] = zz;
698 index += 3;
699 atom.moveTo(xx, yy, zz);
700 }
701 } else {
702
703 molecularAssembly.moveToFractionalCoordinates();
704 index = 0;
705 for (int i = 0; i < nActive; i++) {
706 Atom atom = activeAtoms[i];
707 double xx = atom.getX();
708 double yy = atom.getY();
709 double zz = atom.getZ();
710 xyz[index] = xx;
711 xyz[index + 1] = yy;
712 xyz[index + 2] = zz;
713 index += 3;
714 }
715 }
716 }
717
718 }