n2p2 - A neural network potential package
Loading...
Searching...
No Matches
fix_hdnnp.h
Go to the documentation of this file.
1/* -*- c++ -*- ----------------------------------------------------------
2 LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
3 http://lammps.sandia.gov, Sandia National Laboratories
4 Steve Plimpton, sjplimp@sandia.gov
5
6 Copyright (2003) Sandia Corporation. Under the terms of Contract
7 DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
8 certain rights in this software. This software is distributed under
9 the GNU General Public License.
10
11 See the README file in the top-level LAMMPS directory.
12------------------------------------------------------------------------- */
13#ifdef FIX_CLASS
14
15FixStyle(hdnnp,FixHDNNP)
16
17#else
18
19#ifndef LMP_FIX_HDNNP_H
20#define LMP_FIX_HDNNP_H
21
22#include <gsl/gsl_multimin.h>
23#include "fix.h"
24#include "InterfaceLammps.h"
25
26
27namespace LAMMPS_NS {
28
29 class FixHDNNP : public Fix {
30 friend class PairHDNNP4G;
31 public:
32 FixHDNNP(class LAMMPS *, int, char **);
33 ~FixHDNNP();
34
35 int setmask();
36 virtual void post_constructor();
37 virtual void init();
38 void init_list(int,class NeighList *);
39 void setup_pre_force(int);
40 virtual void pre_force(int);
41
42 void min_setup_pre_force(int);
43 void min_pre_force(int);
44
45 protected:
46
47 class PairHDNNP4G *hdnnp; // interface to HDNNP pair_style
48 class KSpaceHDNNP *kspacehdnnp; // interface to HDNNP kspace_style
49 class NeighList *list;
51
53 int kspaceflag; // 0:Ewald Sum, 1:PPPM
54 int ngroup;
55
56 bool periodic; // true if periodic
57 double qRef; // total reference charge of the system
58 double *Q;
59
60 virtual void pertype_parameters(char*);
61 void isPeriodic(); // check for periodicity
62 void calculate_electronegativities(); // calculates electronegatives via running first set of NNs
63 void process_first_network(); // interfaces to n2p2 and runs first NN
64 void allocate_QEq(); // allocates QEq arrays
65 void deallocate_QEq(); // deallocates QEq arrays
66
68
69 gsl_multimin_function_fdf QEq_minimizer; // minimizer object
70 const gsl_multimin_fdfminimizer_type *T;
71 gsl_multimin_fdfminimizer *s;
72
73 double QEq_f(const gsl_vector*); // f : QEq energy as a function of atomic charges
74 void QEq_df(const gsl_vector*, gsl_vector*); // df : Gradient of QEq energy with respect to atomic charges
75 void QEq_fdf(const gsl_vector*, double*, gsl_vector*); // f * df
76
77 static double QEq_f_wrap(const gsl_vector*, void*); // wrapper of f
78 static void QEq_df_wrap(const gsl_vector*, void*, gsl_vector*); // wrapper of df
79 static void QEq_fdf_wrap(const gsl_vector*, void*, double*, gsl_vector*); // wrapper of f * df
80
81 void calculate_QEqCharges(); // main function where minimization happens
82
83 void calculate_erfc_terms(); // loops over neighbors of local atoms and calculates erfc terms
84 // these consume the most time. Storing them speeds up calculations
85
88 double *qall,*qall_loc;
89 double *dEdQ_all,*dEdQ_loc; // gradient of the charge equilibration energy wrt charges
90 double *xx,*xy,*xz; // global positions (here only for nonperiodic case)
91 double *xx_loc,*xy_loc,*xz_loc; // sparse local positions (here only for nonperiodic case)
92
93
94 /* Matrix Approach (DEPRECATED and to be deleted)
95 int nevery,hdnnpflag;
96 int n, N, m_fill;
97 int n_cap, m_cap;
98 int pack_flag;
99 bigint ngroup;
100 int nprev;
101 double **Q_hist;
102 double *p, *q, *r, *d;
103 virtual int pack_forward_comm(int, int *, double *, int, int *);
104 virtual void unpack_forward_comm(int, int, double *);
105 virtual int pack_reverse_comm(int, int, double *);
106 virtual void unpack_reverse_comm(int, int *, double *);
107 virtual double memory_usage();
108 virtual void grow_arrays(int);
109 virtual void copy_arrays(int, int, int);
110 virtual int pack_exchange(int, double *);
111 virtual int unpack_exchange(int, double *);
112
113 virtual double parallel_norm( double*, int );
114 virtual double parallel_dot( double*, double*, int );
115 virtual double parallel_vector_acc( double*, int );
116
117 virtual void vector_sum(double*,double,double*,double,double*,int);
118 virtual void vector_add(double*, double, double*,int);*/
119
120 };
121
122}
123
124#endif
125#endif
virtual void pertype_parameters(char *)
gsl_multimin_fdfminimizer * s
Definition fix_hdnnp.h:71
class PairHDNNP4G * hdnnp
Definition fix_hdnnp.h:47
void QEq_fdf(const gsl_vector *, double *, gsl_vector *)
void min_setup_pre_force(int)
virtual void pre_force(int)
const gsl_multimin_fdfminimizer_type * T
Definition fix_hdnnp.h:70
friend class PairHDNNP4G
Definition fix_hdnnp.h:30
static double QEq_f_wrap(const gsl_vector *, void *)
static void QEq_df_wrap(const gsl_vector *, void *, gsl_vector *)
void init_list(int, class NeighList *)
FixHDNNP(class LAMMPS *, int, char **)
Definition fix_hdnnp.cpp:52
int * type_all
Global storage.
Definition fix_hdnnp.h:87
class NeighList * list
Definition fix_hdnnp.h:49
gsl_multimin_function_fdf QEq_minimizer
QEq energy minimization via gsl library.
Definition fix_hdnnp.h:69
virtual void post_constructor()
void QEq_df(const gsl_vector *, gsl_vector *)
static void QEq_fdf_wrap(const gsl_vector *, void *, double *, gsl_vector *)
double QEq_f(const gsl_vector *)
void calculate_electronegativities()
class KSpaceHDNNP * kspacehdnnp
Definition fix_hdnnp.h:48
void min_pre_force(int)
virtual void init()
void setup_pre_force(int)