n2p2 - A neural network potential package
Loading...
Searching...
No Matches
LAMMPS_NS::PairHDNNPDevelop Class Reference

#include <pair_hdnnp_develop.h>

Inheritance diagram for LAMMPS_NS::PairHDNNPDevelop:
Collaboration diagram for LAMMPS_NS::PairHDNNPDevelop:

Public Member Functions

 PairHDNNPDevelop (class LAMMPS *)
 
void compute (int, int) override
 
void init_style () override
 
double init_one (int i, int j) override
 
- Public Member Functions inherited from LAMMPS_NS::PairHDNNP
 PairHDNNP (class LAMMPS *)
 
 ~PairHDNNP () override
 
void compute (int, int) override
 
void settings (int, char **) override
 
void coeff (int, char **) override
 
void init_style () override
 
double init_one (int, int) override
 

Protected Member Functions

void transferNeighborList (double const cutoffRadius)
 
void updateNeighborlistCutoff ()
 Update neighborlist if maxCutoffRadiusNeighborList has changed.
 
- Protected Member Functions inherited from LAMMPS_NS::PairHDNNP
virtual void allocate ()
 
void transferNeighborList ()
 
void handleExtrapolationWarnings ()
 

Protected Attributes

double maxCutoffRadiusNeighborList
 Keeps track of the maximum cutoff radius that was used for the neighbor list.
 
- Protected Attributes inherited from LAMMPS_NS::PairHDNNP
bool showew
 
bool resetew
 
int showewsum
 
int maxew
 
long numExtrapolationWarningsTotal
 
long numExtrapolationWarningsSummary
 
double cflength
 
double cfenergy
 
double maxCutoffRadius
 
char * directory
 
std::string emap
 
nnp::InterfaceLammpsinterface
 

Detailed Description

Definition at line 30 of file pair_hdnnp_develop.h.

Constructor & Destructor Documentation

◆ PairHDNNPDevelop()

PairHDNNPDevelop::PairHDNNPDevelop ( class LAMMPS * lmp)

Definition at line 44 of file pair_hdnnp_develop.cpp.

44: PairHDNNP(lmp) {}
PairHDNNP(class LAMMPS *)

References LAMMPS_NS::PairHDNNP::PairHDNNP().

Here is the call graph for this function:

Member Function Documentation

◆ compute()

void PairHDNNPDevelop::compute ( int eflag,
int vflag )
override

Definition at line 48 of file pair_hdnnp_develop.cpp.

49{
50 if(eflag || vflag) ev_setup(eflag,vflag);
51 else evflag = vflag_fdotr = eflag_global = eflag_atom = 0;
52
53 // Set number of local atoms and add index and element.
54 interface->setLocalAtoms(atom->nlocal, atom->type);
55 // Transfer tags separately. Interface::setLocalTags is overloaded internally
56 // to work with both -DLAMMPS_SMALLBIG (tagint = int) and -DLAMMPS_BIGBIG
57 // (tagint = int64_t)
58 interface->setLocalTags(atom->tag);
59
60
61 // Also set absolute atom positions.
62 interface->setLocalAtomPositions(atom->x);
63
64 // Set Box vectors if system is periodic in all 3 dims.
65 if(domain->nonperiodic == 0)
66 {
67 interface->setBoxVectors(domain->boxlo,
68 domain->boxhi,
69 domain->xy,
70 domain->xz,
71 domain->yz);
72 }
73
75
76 // Transfer local neighbor list to NNP interface. Has to be called after
77 // setBoxVectors if structure is periodic!
79
80 // Compute symmetry functions, atomic neural networks and add up energy.
81 interface->processDevelop();
82
83 // Do all stuff related to extrapolation warnings.
84 if(showew == true || showewsum > 0 || maxew >= 0) {
86 }
87
88 // Calculate forces of local and ghost atoms.
89 interface->getForcesDevelop(atom->f);
90
91 // Transfer charges LAMMPS. Does nothing if nnpType != 4G.
92 interface->getCharges(atom->q);
93
94 // Add energy contribution to total energy.
95 if (eflag_global)
96 ev_tally(0,0,atom->nlocal,1,interface->getEnergy(),0.0,0.0,0.0,0.0,0.0);
97
98 // Add atomic energy if requested (CAUTION: no physical meaning!).
99 if (eflag_atom)
100 for (int i = 0; i < atom->nlocal; ++i)
101 eatom[i] = interface->getAtomicEnergy(i);
102
103 // If virial needed calculate via F dot r.
104 if (vflag_fdotr) virial_fdotr_compute();
105
106 //interface.writeToFile("md_run.data", true);
107}
double maxCutoffRadiusNeighborList
Keeps track of the maximum cutoff radius that was used for the neighbor list.
void updateNeighborlistCutoff()
Update neighborlist if maxCutoffRadiusNeighborList has changed.
void handleExtrapolationWarnings()
nnp::InterfaceLammps * interface
Definition pair_hdnnp.h:65

References LAMMPS_NS::PairHDNNP::handleExtrapolationWarnings(), LAMMPS_NS::PairHDNNP::interface, maxCutoffRadiusNeighborList, LAMMPS_NS::PairHDNNP::maxew, LAMMPS_NS::PairHDNNP::showew, LAMMPS_NS::PairHDNNP::showewsum, LAMMPS_NS::PairHDNNP::transferNeighborList(), and updateNeighborlistCutoff().

Here is the call graph for this function:

◆ init_style()

void PairHDNNPDevelop::init_style ( )
override

Definition at line 113 of file pair_hdnnp_develop.cpp.

114{
115 if (comm->nprocs > 1) {
116 throw runtime_error("ERROR: Pair style \"nnp/develop\" can only be used "
117 "with a single MPI task.\n");
118 }
119
120 if(domain->dimension != 3)
121 throw runtime_error("ERROR: Only 3d systems can be used!");
122
123 if (!(domain->xperiodic == domain->yperiodic
124 && domain->yperiodic == domain->zperiodic))
125 throw runtime_error("ERROR: System must be either aperiodic or periodic "
126 "in all 3 dimmensions!");
127
128 // Required for charge equilibration scheme.
129 if (atom->map_style == Atom::MAP_NONE)
130 throw runtime_error("ERROR: pair style requires atom map yes");
131
133
134
136
137 interface->setGlobalStructureStatus(true);
138}
void init_style() override

References LAMMPS_NS::PairHDNNP::init_style(), LAMMPS_NS::PairHDNNP::interface, LAMMPS_NS::PairHDNNP::maxCutoffRadius, and maxCutoffRadiusNeighborList.

Here is the call graph for this function:

◆ init_one()

double PairHDNNPDevelop::init_one ( int i,
int j )
override

Definition at line 141 of file pair_hdnnp_develop.cpp.

142{
143 //cutsq[i][j] = cutsq[j][i] = pow(maxCutoffRadiusNeighborList,2);
145}

References maxCutoffRadiusNeighborList.

◆ transferNeighborList()

void PairHDNNPDevelop::transferNeighborList ( double const cutoffRadius)
protected

Definition at line 148 of file pair_hdnnp_develop.cpp.

149{
150 // Transfer neighbor list to NNP.
151 double rc2 = cutoffRadius * cutoffRadius;
152 interface->allocateNeighborlists(list->numneigh);
153#ifdef _OPENMP
154 #pragma omp parallel for
155#endif
156 for (int ii = 0; ii < list->inum; ++ii) {
157 int i = list->ilist[ii];
158 for (int jj = 0; jj < list->numneigh[i]; ++jj) {
159 int j = list->firstneigh[i][jj];
160 j &= NEIGHMASK;
161 double dx = atom->x[i][0] - atom->x[j][0];
162 double dy = atom->x[i][1] - atom->x[j][1];
163 double dz = atom->x[i][2] - atom->x[j][2];
164 double d2 = dx * dx + dy * dy + dz * dz;
165 if (d2 <= rc2) {
166 if (!interface->getGlobalStructureStatus())
167 // atom->tag[j] will be implicitly converted to int64_t internally.
168 interface->addNeighbor(i,j,atom->tag[j],atom->type[j],dx,dy,dz,d2);
169 else
170 interface->addNeighbor(i,j,atom->map(atom->tag[j]),atom->type[j],dx,dy,dz,d2);
171 }
172 }
173 }
174 interface->finalizeNeighborList();
175}

References LAMMPS_NS::PairHDNNP::interface.

◆ updateNeighborlistCutoff()

void PairHDNNPDevelop::updateNeighborlistCutoff ( )
protected

Update neighborlist if maxCutoffRadiusNeighborList has changed.

Definition at line 178 of file pair_hdnnp_develop.cpp.

179{
180 double maxCutoffRadiusOverall = interface->getMaxCutoffRadiusOverall();
181 if(maxCutoffRadiusOverall > maxCutoffRadiusNeighborList)
182 {
183 // TODO: Increase slightly to compensate for rounding errors?
184 utils::logmesg(lmp, fmt::format("WARNING: Cutoff too small, need at "
185 "least: {}\n", maxCutoffRadiusOverall));
186
187 //maxCutoffRadiusNeighborList = maxCutoffRadiusOverall;
188 //reinit();
189
190 //mixed_flag = 1;
191 //for (int i = 1; i <= atom->ntypes; i++)
192 //{
193 // for (int j = i; j <= atom->ntypes; j++)
194 // {
195 // if ((i != j) && setflag[i][j]) mixed_flag = 0;
196 // cutsq[i][j] = cutsq[j][i] = pow(maxCutoffRadiusNeighborList, 2);
197 // }
198 //}
199 //cutforce = MAX(cutforce, maxCutoffRadiusNeighborList);
200
201 //neighbor->init();
202 //neighbor->build(0);
203 }
204}

References LAMMPS_NS::PairHDNNP::interface, and maxCutoffRadiusNeighborList.

Referenced by compute().

Here is the caller graph for this function:

Member Data Documentation

◆ maxCutoffRadiusNeighborList

double LAMMPS_NS::PairHDNNPDevelop::maxCutoffRadiusNeighborList
protected

Keeps track of the maximum cutoff radius that was used for the neighbor list.

Definition at line 41 of file pair_hdnnp_develop.h.

Referenced by compute(), init_one(), init_style(), and updateNeighborlistCutoff().


The documentation for this class was generated from the following files: