runtime schedule adding and injection of c atoms
[physik/posic.git] / moldyn.h
index b140372..51ad482 100644 (file)
--- a/moldyn.h
+++ b/moldyn.h
@@ -77,6 +77,7 @@ typedef struct s_moldyn_schedule {
 /* moldyn main structure */
 typedef struct s_moldyn {
        int count;              /* total amount of atoms */
+       double mass;            /* total system mass */
        t_atom *atom;           /* pointer to the atoms */
 
        t_3dvec dim;            /* dimensions of the simulation volume */
@@ -85,11 +86,15 @@ typedef struct s_moldyn {
        /* potential force function and parameter pointers */
        int (*func1b)(struct s_moldyn *moldyn,t_atom *ai);
        int (*func2b)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc);
-       int (*func2b_post)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc);
-       int (*func3b)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,t_atom *ak,
-                     u8 bck);
+       int (*func3b_j1)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc);
+       int (*func3b_j2)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc);
+       int (*func3b_j3)(struct s_moldyn *moldyn,t_atom *ai,t_atom *aj,u8 bc);
+       int (*func3b_k1)(struct s_moldyn *moldyn,
+                        t_atom *ai,t_atom *aj,t_atom *ak,u8 bck);
+       int (*func3b_k2)(struct s_moldyn *moldyn,
+                        t_atom *ai,t_atom *aj,t_atom *ak,u8 bck);
        void *pot_params;
-       //int (*potential_force_function)(struct s_moldyn *moldyn);
+       unsigned char run3bp;
 
        double cutoff;          /* cutoff radius */
        double cutoff_square;   /* square of the cutoff radius */
@@ -97,11 +102,31 @@ typedef struct s_moldyn {
 
        t_linkcell lc;          /* linked cell list interface */
 
+       int avg_skip;           /* amount of steps without average calc */
+
        double t_ref;           /* reference temperature */
        double t;               /* actual temperature */
+       double t_sum;           /* sum over all t */
+       double t_avg;           /* average value of t */
+
+       t_virial gvir;          /* global virial (absolute coordinates) */
+       double gv;
+       double gv_sum;
+       double gv_avg;
+
+       double gp;              /* pressure computed from global virial */
+       double gp_sum;          /* sum over all gp */
+       double gp_avg;          /* average value of gp */
+
+       double virial;          /* actual virial */
+       double virial_sum;      /* sum over all calculated virials */
+       double virial_avg;      /* average of virial */
 
        double p_ref;           /* reference pressure */
        double p;               /* actual pressure (computed by virial) */
+       double p_sum;           /* sum over all p */
+       double p_avg;           /* average value of p */
+
        t_3dvec tp;             /* thermodynamic pressure dU/dV */
        double dv;              /* dV for thermodynamic pressure calc */
 
@@ -121,11 +146,28 @@ typedef struct s_moldyn {
        double tau;             /* delta t */
        double time;            /* absolute time */
        double tau_square;      /* delta t squared */
-       double elapsed;         /* total elapsed time */
+       int total_steps;        /* total steps */
 
+       /* energy */
        double energy;          /* potential energy */
        double ekin;            /* kinetic energy */
 
+       /* energy averages & fluctuations */
+       double k_sum;           /* sum of kinetic energy */
+       double v_sum;           /* sum of potential energy */
+       double k_avg;           /* average of kinetic energy */
+       double v_avg;           /* average of potential energy */
+       double k2_sum;          /* sum of kinetic energy squared */
+       double v2_sum;          /* sum of potential energy squared */
+       double k2_avg;          /* average of kinetic energy squared */
+       double v2_avg;          /* average of potential energy squared */
+       double dk2_avg;         /* mean square kinetic energy fluctuations */
+       double dv2_avg;         /* mean square potential energy fluctuations */
+       
+       /* response functions */
+       double c_v_nve;         /* constant volume heat capacity (nve) */
+       double c_v_nvt;         /* constant volume heat capacity (nvt) */
+
        char vlsdir[128];       /* visualization/log/save directory */
        t_visual vis;           /* visualization interface structure */
        u8 vlsprop;             /* log/vis/save properties */
@@ -133,12 +175,18 @@ typedef struct s_moldyn {
        int efd;                /* fd for energy log */
        unsigned int mwrite;    /* how often to log momentum */
        int mfd;                /* fd for momentum log */
+       unsigned int pwrite;    /* how often to log pressure */
+       int pfd;                /* fd for pressure log */
+       unsigned int twrite;    /* how often to log temperature */
+       int tfd;                /* fd for temperature log */
        unsigned int vwrite;    /* how often to visualize atom information */
        unsigned int swrite;    /* how often to create a save file */
        int rfd;                /* report file descriptor */
        char rtitle[64];        /* report title */
        char rauthor[64];       /* report author */
-       int pfd;                /* gnuplot script file descriptor */
+       int epfd;               /* energy gnuplot script file descriptor */
+       int ppfd;               /* pressure gnuplot script file descriptor */
+       int tpfd;               /* temperature gnuplot script file descriptor */
 
        u8 status;              /* general moldyn properties */
 
@@ -169,7 +217,7 @@ typedef struct s_moldyn {
 #define P_SCALE_DIRECT                 0x08    /* direct p control */
 
 /*
- * default values
+ * default values & units
  *
  * - length unit: 1 A (1 A = 1e-10 m)
  * - time unit: 1 fs (1 fs = 1e-15 s)
@@ -185,7 +233,11 @@ typedef struct s_moldyn {
 #define KILOGRAM                       (1.0/AMU)               /* amu */
 #define NEWTON (METER*KILOGRAM/(SECOND*SECOND))        /* A amu / fs^2 */
 #define PASCAL (NEWTON/(METER*METER))                  /* N / A^2 */
-#define ATM    ((1.0133e5*PASCAL))                     /* N / A^2 */
+#define BAR    ((1.0e5*PASCAL))                        /* N / A^2 */
+#define K_BOLTZMANN    (1.380650524e-23*METER*NEWTON)  /* NA/K */
+#define K_B2           (K_BOLTZMANN*K_BOLTZMANN)       /* (NA)^2/K^2 */
+#define EV             (1.6021765314e-19*METER*NEWTON) /* NA */
+#define JOULE          (NEWTON*METER)                  /* NA */
 
 #define MOLDYN_TEMP                    273.0
 #define MOLDYN_TAU                     1.0
@@ -201,9 +253,11 @@ typedef struct s_moldyn {
 
 #define LOG_TOTAL_ENERGY               0x01
 #define LOG_TOTAL_MOMENTUM             0x02
-#define SAVE_STEP                      0x04
-#define VISUAL_STEP                    0x08
-#define CREATE_REPORT                  0x10
+#define LOG_PRESSURE                   0x04
+#define LOG_TEMPERATURE                        0x08
+#define SAVE_STEP                      0x10
+#define VISUAL_STEP                    0x20
+#define CREATE_REPORT                  0x40
 
 #define TRUE                           1
 #define FALSE                          0
@@ -211,48 +265,50 @@ typedef struct s_moldyn {
 #define VERBOSE                                1
 #define QUIET                          0
 
+#define SCALE_UP                       'u'
+#define SCALE_DOWN                     'd'
+#define SCALE_DIRECT                   'D'
+
 /*
- *
- * phsical values / constants
- *
+ * potential related phsical values / constants
  *
  */
 
 #define ONE_THIRD              (1.0/3.0)
 
-#define K_BOLTZMANN            (1.380650524e-23*METER*NEWTON)  /* NA/K */
-#define EV                     (1.6021765314e-19*METER*NEWTON) /* NA */
-
 #define C                      0x06
+#define LC_C                   3.567                           /* A */
 #define M_C                    12.011                          /* amu */
 
 #define SI                     0x0e
-#define LC_SI                  (0.543105e-9*METER)             /* A */
+#define LC_SI                  5.43105                         /* A */
 #define M_SI                   28.08553                        /* amu */
 
+#define LC_3C_SIC              4.3596                          /* A */
+
 #define LJ_SIGMA_SI            ((0.25*sqrt(3.0)*LC_SI)/1.122462)       /* A */
 //#define LJ_SIGMA_SI          (LC_SI/1.122462)                        /* A */
-//#define LJ_SIGMA_SI          (0.5*sqrt(2.0)*LC_SI/1.122462)                  /* A */
+//#define LJ_SIGMA_SI          (0.5*sqrt(2.0)*LC_SI/1.122462)          /* A */
 #define LJ_EPSILON_SI          (2.1678*EV)                             /* NA */
 
-#define TM_R_SI                        (2.7e-10*METER)                 /* A */
-#define TM_S_SI                        (3.0e-10*METER)                 /* A */
+#define TM_R_SI                        2.7                             /* A */
+#define TM_S_SI                        3.0                             /* A */
 #define TM_A_SI                        (1830.8*EV)                     /* NA */
 #define TM_B_SI                        (471.18*EV)                     /* NA */
-#define TM_LAMBDA_SI           (2.4799e10/METER)               /* 1/A */
-#define TM_MU_SI               (1.7322e10/METER)               /* 1/A */
+#define TM_LAMBDA_SI           2.4799                          /* 1/A */
+#define TM_MU_SI               1.7322                          /* 1/A */
 #define TM_BETA_SI             1.1000e-6
 #define TM_N_SI                        0.78734
 #define TM_C_SI                        1.0039e5
 #define TM_D_SI                        16.217
 #define TM_H_SI                        -0.59825
 
-#define TM_R_C                 (1.8e-10*METER)                 /* A */
-#define TM_S_C                 (2.1e-10*METER)                 /* A */
+#define TM_R_C                 1.8                             /* A */
+#define TM_S_C                 2.1                             /* A */
 #define TM_A_C                 (1393.6*EV)                     /* NA */
 #define TM_B_C                 (346.7*EV)                      /* NA */
-#define TM_LAMBDA_C            (3.4879e10/METER)               /* 1/A */
-#define TM_MU_C                        (2.2119e10/METER)               /* 1/A */
+#define TM_LAMBDA_C            3.4879                          /* 1/A */
+#define TM_MU_C                        2.2119                          /* 1/A */
 #define TM_BETA_C              1.5724e-7
 #define TM_N_C                 0.72751
 #define TM_C_C                 3.8049e4
@@ -261,8 +317,53 @@ typedef struct s_moldyn {
 
 #define TM_CHI_SIC             0.9776
 
+#define TM_LC_SIC              4.32                            /* A */
+
+#define ALBE_R_SI              (2.82-0.14)
+#define ALBE_S_SI              (2.82+0.14)
+#define ALBE_A_SI              (3.24*EV/0.842)
+#define ALBE_B_SI              (1.842*3.24*EV/0.842)
+#define ALBE_R0_SI             2.232
+#define ALBE_LAMBDA_SI         (1.4761*sqrt(2.0*1.842))
+#define ALBE_MU_SI             (1.4761*sqrt(2.0/1.842))
+#define ALBE_GAMMA_SI          0.114354
+#define ALBE_C_SI              2.00494
+#define ALBE_D_SI              0.81472
+#define ALBE_H_SI              0.259
+
+#define ALBE_LC_SI             5.429
+
+#define ALBE_R_C               (2.00-0.15)
+#define ALBE_S_C               (2.00+0.15)
+#define ALBE_A_C               (6.00*EV/1.167)
+#define ALBE_B_C               (2.167*6.00*EV/1.167)
+#define ALBE_R0_C              1.4276
+#define ALBE_LAMBDA_C          (2.0099*sqrt(2.0*2.167))
+#define ALBE_MU_C              (2.0099*sqrt(2.0/2.167))
+#define ALBE_GAMMA_C           0.11233
+#define ALBE_C_C               181.910
+#define ALBE_D_C               6.28433
+#define ALBE_H_C               0.5556
+
+#define ALBE_LC_C              3.566
+
+#define ALBE_R_SIC             (2.40-0.20)
+#define ALBE_S_SIC             (2.40+0.10)
+#define ALBE_A_SIC             (4.36*EV/0.847)
+#define ALBE_B_SIC             (1.847*4.36*EV/0.847)
+#define ALBE_R0_SIC            1.79
+#define ALBE_LAMBDA_SIC                (1.6991*sqrt(2.0*1.847))
+#define ALBE_MU_SIC            (1.6991*sqrt(2.0/1.847))
+#define ALBE_GAMMA_SIC         0.011877
+#define ALBE_C_SIC             273987
+#define ALBE_D_SIC             180.314
+#define ALBE_H_SIC             0.68
+
+#define ALBE_LC_SIC            4.359
+
+
 /*
- * lattice constants
+ * lattice types
  */
 
 #define CUBIC                  0x01
@@ -276,10 +377,9 @@ typedef struct s_moldyn {
  *
  */
 
-typedef int (*pf_func1b)(t_moldyn *,t_atom *ai);
-typedef int (*pf_func2b)(t_moldyn *,t_atom *,t_atom *,u8 bc);
-typedef int (*pf_func2b_post)(t_moldyn *,t_atom *,t_atom *,u8 bc);
-typedef int (*pf_func3b)(t_moldyn *,t_atom *,t_atom *,t_atom *,u8 bc);
+typedef int (*pf_func1b)(t_moldyn *,t_atom *);
+typedef int (*pf_func2b)(t_moldyn *,t_atom *,t_atom *,u8);
+typedef int (*pf_func3b)(t_moldyn *,t_atom *,t_atom *,t_atom *,u8);
 
 int moldyn_init(t_moldyn *moldyn,int argc,char **argv);
 int moldyn_shutdown(t_moldyn *moldyn);
@@ -294,17 +394,22 @@ int set_nn_dist(t_moldyn *moldyn,double dist);
 int set_pbc(t_moldyn *moldyn,u8 x,u8 y,u8 z);
 int set_potential1b(t_moldyn *moldyn,pf_func1b func);
 int set_potential2b(t_moldyn *moldyn,pf_func2b func);
-int set_potential2b_post(t_moldyn *moldyn,pf_func2b_post func);
-int set_potential3b(t_moldyn *moldyn,pf_func3b func);
+int set_potential3b_j1(t_moldyn *moldyn,pf_func2b func);
+int set_potential3b_j2(t_moldyn *moldyn,pf_func2b func);
+int set_potential3b_j3(t_moldyn *moldyn,pf_func2b func);
+int set_potential3b_k1(t_moldyn *moldyn,pf_func3b func);
+int set_potential3b_k2(t_moldyn *moldyn,pf_func3b func);
 int set_potential_params(t_moldyn *moldyn,void *params);
 
+int set_avg_skip(t_moldyn *moldyn,int skip);
+
 int moldyn_set_log_dir(t_moldyn *moldyn,char *dir);
 int moldyn_set_report(t_moldyn *moldyn,char *author,char *title);
 int moldyn_set_log(t_moldyn *moldyn,u8 type,int timer);
 int moldyn_log_shutdown(t_moldyn *moldyn);
 
 int create_lattice(t_moldyn *moldyn,u8 type,double lc,int element,double mass,
-                   u8 attr,u8 brand,int a,int b,int c);
+                   u8 attr,u8 brand,int a,int b,int c,t_3dvec *origin);
 int cubic_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin);
 int fcc_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin);
 int diamond_init(int a,int b,int c,double lc,t_atom *atom,t_3dvec *origin);
@@ -313,18 +418,21 @@ int add_atom(t_moldyn *moldyn,int element,double mass,u8 brand,u8 attr,
 int destroy_atoms(t_moldyn *moldyn);
 
 int thermal_init(t_moldyn *moldyn,u8 equi_init);
+double total_mass_calc(t_moldyn *moldyn);
 double temperature_calc(t_moldyn *moldyn);
 double get_temperature(t_moldyn *moldyn);
 int scale_velocity(t_moldyn *moldyn,u8 equi_init);
+double virial_sum(t_moldyn *moldyn);
 double pressure_calc(t_moldyn *moldyn);
+int energy_fluctuation_calc(t_moldyn *moldyn);
+int get_heat_capacity(t_moldyn *moldyn);
 double thermodynamic_pressure_calc(t_moldyn *moldyn);
 double get_pressure(t_moldyn *moldyn);
 int scale_volume(t_moldyn *moldyn);
-int scale_dim(t_moldyn *moldyn,double scale,u8 x,u8 y,u8 z);
-int scale_atoms(t_moldyn *moldyn,double scale,u8 x,u8 y,u8 z);
+int scale_dim(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z);
+int scale_atoms(t_moldyn *moldyn,u8 dir,double scale,u8 x,u8 y,u8 z);
 
-double get_e_kin(t_moldyn *moldyn);
-double update_e_kin(t_moldyn *moldyn);
+double e_kin_calc(t_moldyn *moldyn);
 double get_total_energy(t_moldyn *moldyn);
 t_3dvec get_total_p(t_moldyn *moldyn);
 
@@ -353,4 +461,6 @@ int check_per_bound(t_moldyn *moldyn,t_3dvec *a);
 
 int moldyn_bc_check(t_moldyn *moldyn);
 
+int get_line(int fd,char *line,int max);
+
 #endif