1 /* ----------------------------------------------------------------------
2     This is the
3 
4     ██╗     ██╗ ██████╗  ██████╗  ██████╗ ██╗  ██╗████████╗███████╗
5     ██║     ██║██╔════╝ ██╔════╝ ██╔════╝ ██║  ██║╚══██╔══╝██╔════╝
6     ██║     ██║██║  ███╗██║  ███╗██║  ███╗███████║   ██║   ███████╗
7     ██║     ██║██║   ██║██║   ██║██║   ██║██╔══██║   ██║   ╚════██║
8     ███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║  ██║   ██║   ███████║
9     ╚══════╝╚═╝ ╚═════╝  ╚═════╝  ╚═════╝ ╚═╝  ╚═╝   ╚═╝   ╚══════╝®
10 
11     DEM simulation engine, released by
12     DCS Computing Gmbh, Linz, Austria
13     http://www.dcs-computing.com, office@dcs-computing.com
14 
15     LIGGGHTS® is part of CFDEM®project:
16     http://www.liggghts.com | http://www.cfdem.com
17 
18     Core developer and main author:
19     Christoph Kloss, christoph.kloss@dcs-computing.com
20 
21     LIGGGHTS® is open-source, distributed under the terms of the GNU Public
22     License, version 2 or later. It is distributed in the hope that it will
23     be useful, but WITHOUT ANY WARRANTY; without even the implied warranty
24     of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have
25     received a copy of the GNU General Public License along with LIGGGHTS®.
26     If not, see http://www.gnu.org/licenses . See also top-level README
27     and LICENSE files.
28 
29     LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH,
30     the producer of the LIGGGHTS® software and the CFDEM®coupling software
31     See http://www.cfdem.com/terms-trademark-policy for details.
32 
33 -------------------------------------------------------------------------
34     Contributing author and copyright for this file:
35     Christoph Kloss (JKU Linz, DCS Computing GmbH, Linz)
36     Richard Berger (JKU Linz)
37 
38     Copyright 2012-     DCS Computing GmbH, Linz
39     Copyright 2009-2015 JKU Linz
40 ------------------------------------------------------------------------- */
41 
42 #ifdef FIX_CLASS
43 
44 #else
45 
46 #ifndef LMP_FIX_INSERT_H
47 #define LMP_FIX_INSERT_H
48 
49 #include "fix.h"
50 #include "bounding_box.h"
51 #include "region_neighbor_list.h"
52 #include "fix_particledistribution_discrete.h"
53 #include "update.h"
54 
55 namespace LAMMPS_NS {
56 
57 class FixInsert : public Fix {
58  public:
59   FixInsert(class LAMMPS *, int, char **);
60   ~FixInsert();
61 
62   virtual int setmask();
63   virtual void init();
64   void reset_timestep(bigint newstep,bigint oldstep);
setup_pre_exchange()65   virtual void setup_pre_exchange() {}
66   void setup(int vflag);
67   virtual double extend_cut_ghost();
68   void pre_exchange();
end_of_step()69   virtual void end_of_step() {}
70 
71   void write_restart(FILE *);
72   virtual void restart(char *);
73 
reset_timestep(bigint)74   virtual void reset_timestep(bigint) {}
75 
76   double compute_vector(int index);
77 
78   virtual double min_rad(int);
79   virtual double max_rad(int);
80   virtual double max_r_bound();
81   int min_type();
82   int max_type();
83 
ins_every()84   int ins_every()
85   { return insert_every; }
86 
get_distribution()87   FixParticledistributionDiscrete * get_distribution()
88   { return fix_distribution; }
89 
has_inserted()90   bool has_inserted() const
91   { return most_recent_ins_step == update->ntimestep; }
92 
93  protected:
94 
95   int iarg;
96 
97   /*---TARGET QUANTITIES---how much will be inserted*/
98 
99   // flag if total # of particle to insert by this command exists
100   int ninsert_exists;
101 
102   // insertion target - can bei either number or mass
103   // in case of multi-sphere, ninsert is number of bodies to insert
104   int ninsert;
105   double massinsert;
106 
107   // insertion rate - can bei either number rate or mass rate
108   // in particle per time and mass per time units - both can be 0
109   double nflowrate;
110   double massflowrate;
111 
112   // how much we have inserted to far
113   int ninserted;
114   double massinserted;
115 
116   // max allocation for insertion
117   int ninsert_this_max_local;
118 
119   // first, most recent, final insertion step
120   int first_ins_step,most_recent_ins_step, final_ins_step;
121 
122   /*---INSERTION QUANTITIES---what, where and how exactly will we insert*/
123 
124   //particle distribution
125   class FixParticledistributionDiscrete *fix_distribution;
126 
127   // insert 'ninsert_per' particles every 'insert_every' steps
128   // 'ninsert_per' is the default, actual # of inserted particles
129   // to insert is defined in calc_ninsert_this()
130   int insert_every;
131   double ninsert_per;
132 
133   // determines how particle distributions are interpreted
134   // total number of particles inserted for each vary with calc_ninsert_this() result
135   // if exact_number = 0, total # inserted particles will deviate from calc_ninsert_this() result
136   //                      distribution will be fulfilled statistically
137   // if exact_number = 1, total # inserted particles will match result from calc_ninsert_this()
138   int exact_number;
139 
140   // max and min type to be inserted
141   int type_max,type_min;
142 
143   // minmum/maximum radius to be inserted
144   double minrad, maxrad;
145 
146   // flag if overlap is checked upon insertion (via all-to-all comm)
147   int check_ol_flag;
148 #ifdef SUPERQUADRIC_ACTIVE_FLAG
149   int check_obb_flag;
150 #endif
151 
152   // if flag is 1, particles are generated to be in the region as a whole
153   // if flag is 0, particles centers are in region
154   int all_in_flag;
155 
156   int maxattempt;
157 
158   // positions generated, and for overlap check
159   RegionNeighborList<interpolate_no> &neighList;
160 
161   // velocity and ang vel distribution
162   // currently constant for omega - could also be a distribution
163   int    v_randomSetting;
164   double v_insert[3];
165   double v_insertFluct[3];
166   double omega_insert[3];
167   double quat_insert[4];
168   bool quat_random_;
169 
170   /*---FURTHER THINGS THAT WE NEED---*/
171 
172   // random generator
173   class RanPark *random;
174   int seed;
175 
176   // all2all comm
177   int me,nprocs;
178   int *recvcounts,*displs;
179 
180   // determine if print stats
181   int print_stats_start_flag;
182   int print_stats_during_flag;
183 
184   // warn if box extent too small for insertion
185   bool warn_boxentent;
186 
187   // compress atom tags upon insertion
188   bool compress_flag;
189 
190   class FixMultisphere *fix_multisphere;
191   class Multisphere *multisphere;
192 
193   /*---FUNCTION MEMBERS---*/
194 
195   virtual void print_stats_start();
196   virtual void print_stats_during(int,double);
197 
198   virtual void init_defaults();
199   virtual void sanity_check();
200   virtual void calc_insertion_properties() = 0;
201 
pre_insert()202   virtual bool pre_insert() { return true; }
203   virtual int calc_ninsert_this();
204   virtual int load_xnear(int);
205   virtual int is_nearby(int) = 0;
206   virtual BoundingBox getBoundingBox() = 0;
207 
init_list(const int ninsert_this_local)208   virtual void init_list(const int ninsert_this_local)
209   { fix_distribution->random_init_list(ninsert_this_local); }
generate_list(const int ninsert_this_local,const int groupbit,const int exact_number)210   virtual int generate_list(const int ninsert_this_local, const int groupbit, const int exact_number)
211   { return fix_distribution->randomize_list(ninsert_this_local, groupbit, exact_number); }
212   virtual void x_v_omega(int,int&,int&,double&) = 0;
213   virtual double insertion_fraction() = 0;
214 
finalize_insertion(int)215   virtual void finalize_insertion(int){};
216 
has_set_property()217   bool has_set_property() const
218   { return property_name != NULL && fix_property != NULL; }
219 
220  protected:
221   void generate_random_velocity(double * velocity);
222 
223  private:
224 
225   char *property_name;
226   class FixPropertyAtom *fix_property;
227   double fix_property_value;
228 
229   bool setup_flag;
230 
231   class Irregular *irregular;
232 
233   virtual int distribute_ninsert_this(int);
234 };
235 
236 }
237 
238 #endif
239 #endif
240