NAMD
GridForceGrid.C
Go to the documentation of this file.
1 
7 #include <iostream>
8 #include <typeinfo>
9 
10 #include "GridForceGrid.h"
11 #include "Vector.h"
12 #include "Node.h"
13 #include "SimParameters.h"
14 #include "Molecule.h"
15 #include "InfoStream.h"
16 #include "common.h"
17 #include "ComputeGridForce.h"
18 
19 #include "MGridforceParams.h"
20 
21 #define MIN_DEBUG_LEVEL 4
22 //#define DEBUGM
23 #include "Debug.h"
24 
25 #include "GridForceGrid.inl"
26 
27 
28 /*****************/
29 /* GRIDFORCEGRID */
30 /*****************/
31 
32 // Class methods
33 
34 GridforceGrid * GridforceGrid::new_grid(int gridnum, char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams)
35 {
36  GridforceGrid *grid = NULL;
37  if (mgridParams->gridforceLite) {
38  grid = new GridforceLiteGrid(gridnum);
39  } else {
40  grid = new GridforceFullMainGrid(gridnum);
41  }
42 
43  grid->initialize(potfilename, simParams, mgridParams);
44 
45  return grid;
46 }
47 
49 
51 {
52  // Abstract interface for packing a grid into a message. This
53  // could easily be a non-static function as it was before, but is
54  // static to make it similar to unpack_grid below, which does need
55  // to be static since it creates a GridforceGrid object.
56  msg->put(grid->get_grid_type());
57  grid->pack(msg);
58 }
59 
61 {
62  // Abstract interface for unpacking a grid from a message.
63  GridforceGrid *grid = NULL;
64  int type;
65 
66  msg->get(type);
67 
68  switch (type) {
70  grid = new GridforceFullMainGrid(gridnum);
71  break;
73  grid = new GridforceLiteGrid(gridnum);
74  break;
75  default:
76  NAMD_bug("GridforceGrid::unpack_grid called with unknown grid type!");
77  }
78 
79  grid->unpack(msg);
80 
81  return grid;
82 }
83 
85 {
86  // Idea: Go through each grid corner and wrap it to the grid center;
87  // if the position moves, then the grid is too big and we return false
88  DebugM(4, "Checking whether grid fits periodic cell\n");
89  Position center = get_center();
90  for (int i = 0; i < 8; i++) {
91  Position pos = get_corner(i);
92  Position pos_wrapped = wrap_position(pos, lattice);
93  if ((pos - pos_wrapped).length() > 1.) {
94  DebugM(5, "(" << pos << ") != (" << pos_wrapped << ")\n" << endi);
95  return false;
96  }
97  }
98  return true;
99 }
100 
102 {
103  // idx -> (x,y,z) (cell basis coordinates)
104  // 0 -> (0,0,0)
105  // 1 -> (1,0,0)
106  // 2 -> (0,1,0)
107  // 3 -> (1,1,0)
108  // 4 -> (0,0,1)
109  // 5 -> (1,0,1)
110  // 6 -> (0,1,1)
111  // 7 -> (1,1,1)
112  Position pos;
113  if (idx >= 8 || idx < 0) {
114  // invalid index
115  pos = Vector(); // default value of Vector() is strange enough to be a decent "undefined" value (-99999, -99999, -999999)
116  } else if (corners[idx] != Vector()) {
117  // use cached value if possible
118  pos = corners[idx];
119  } else {
120  // must calculate
121  Tensor e = get_e();
122  pos = get_origin();
123  if (idx & (1 << 0)) pos += e * Vector(get_k0()-1, 0, 0);
124  if (idx & (1 << 1)) pos += e * Vector(0, get_k1()-1, 0);
125  if (idx & (1 << 2)) pos += e * Vector(0, 0, get_k2()-1);
126  corners[idx] = pos; // cache for future use
127  DebugM(4, "corner " << idx << " = " << pos << "\n" << endi);
128  }
129  return pos;
130 }
131 
132 
133 /*************************/
134 /* GRIDFORCEFULLBASEGRID */
135 /*************************/
136 
138 {
139  cont[0] = cont[1] = cont[2] = FALSE;
140  grid = NULL;
141  numSubgrids = 0;
142  subgrids = NULL;
143 }
144 
146 {
147  delete[] grid;
148  for (int i = 0; i < numSubgrids; i++) {
149  delete subgrids[i];
150  }
151  delete[] subgrids;
152 }
153 
154 
156 {
157  DebugM(2, "Packing message\n" << endi);
158 
159  msg->put(numSubgrids);
160  msg->put(generation);
161 
162  msg->put(3*sizeof(int), (char*)k);
163  msg->put(3*sizeof(int), (char*)k_nopad);
164  msg->put(size);
165  msg->put(size_nopad);
166  msg->put(3*sizeof(long int), (char*)dk);
167  msg->put(3*sizeof(long int), (char*)dk_nopad);
168  msg->put(factor);
169 
170  msg->put(sizeof(Vector), (char*)&origin);
171  msg->put(sizeof(Vector), (char*)&center);
172  msg->put(sizeof(Tensor), (char*)&e);
173  msg->put(sizeof(Tensor), (char*)&inv);
174 
175 // msg->put(3*sizeof(float), (char*)pad_p);
176 // msg->put(3*sizeof(float), (char*)pad_n);
177  msg->put(3*sizeof(Bool), (char*)cont);
178  msg->put(3*sizeof(float), (char*)offset);
179  msg->put(3*sizeof(float), (char*)gap);
180  msg->put(3*sizeof(float), (char*)gapinv);
181  msg->put(sizeof(Vector), (char*)&scale);
182  msg->put(sizeof(Bool), (char*)&checksize);
183 
184  DebugM(2, "Packing grid, size = " << size << "\n" << endi);
185 
186  msg->put(size*sizeof(float), (char*)grid);
187 
188  DebugM(2, "Packing subgrids\n" << endi);
189 
190  for (int i = 0; i < numSubgrids; i++) {
191  subgrids[i]->pack(msg);
192  }
193 }
194 
196 {
197  DebugM(3, "Unpacking message\n" << endi);
198 // iout << iINFO << CkMyPe() << " Unpacking message\n" << endi;
199 
200  delete[] grid;
201  grid = NULL;
202  for (int i = 0; i < numSubgrids; i++) {
203  delete subgrids[i];
204  }
205  numSubgrids = 0;
206  delete[] subgrids;
207  subgrids = NULL;
208 
209  msg->get(numSubgrids);
210  msg->get(generation);
211 
212  DebugM(3, "numSubgrids = " << numSubgrids << "\n");
213  DebugM(3, "generation = " << generation << "\n" << endi);
214 
215  msg->get(3*sizeof(int), (char*)k);
216  msg->get(3*sizeof(int), (char*)k_nopad);
217  msg->get(size);
218  msg->get(size_nopad);
219  msg->get(3*sizeof(long int), (char*)dk);
220  msg->get(3*sizeof(long int), (char*)dk_nopad);
221  msg->get(factor);
222 
223  DebugM(3, "size = " << size << "\n" << endi);
224 
225  msg->get(sizeof(Vector), (char*)&origin);
226  msg->get(sizeof(Vector), (char*)&center);
227  msg->get(sizeof(Tensor), (char*)&e);
228  msg->get(sizeof(Tensor), (char*)&inv);
229 
230 // msg->get(3*sizeof(float), (char*)pad_p);
231 // msg->get(3*sizeof(float), (char*)pad_n);
232  msg->get(3*sizeof(Bool), (char*)cont);
233  msg->get(3*sizeof(float), (char*)offset);
234  msg->get(3*sizeof(float), (char*)gap);
235  msg->get(3*sizeof(float), (char*)gapinv);
236  msg->get(sizeof(Vector), (char*)&scale);
237  msg->get(sizeof(Bool), (char*)&checksize);
238 
239  if (size) {
240  DebugM(3, "allocating grid, size = " << size << "\n" << endi);
241  grid = new float[size];
242  msg->get(size*sizeof(float), (char*)grid);
243  }
244 
245  if (numSubgrids) {
246  DebugM(3, "Creating subgrids array, size " << numSubgrids << "\n" << endi);
248  for (int i = 0; i < numSubgrids; i++) {
249  subgrids[i] = new GridforceFullSubGrid(this);
250  subgrids[i]->unpack(msg);
251  }
252  }
253 }
254 
256 {
257  char line[256];
258  long int poten_offset;
259  do {
260  poten_offset = ftell(poten_fp);
261  fgets(line, 256, poten_fp); // Read comment lines
262  DebugM(4, "Read line: " << line << endi);
263  } while (line[0] == '#');
264  fseek(poten_fp, poten_offset, SEEK_SET);
265 
266  // read grid dimensions
267  fscanf(poten_fp, "object %*d class gridpositions counts %d %d %d\n",
268  &k_nopad[0], &k_nopad[1], &k_nopad[2]);
269  size_nopad = k_nopad[0] * k_nopad[1] * k_nopad[2];
270 
271  // Read origin
272  fscanf(poten_fp, "origin %lf %lf %lf\n",
273  &origin.x, &origin.y, &origin.z);
274 
275  // Read delta (unit vectors)
276  // These are column vectors, so must fill gridfrcE tensor to reflect this
277  fscanf(poten_fp, "delta %lf %lf %lf\n", &e.xx, &e.yx, &e.zx);
278  fscanf(poten_fp, "delta %lf %lf %lf\n", &e.xy, &e.yy, &e.zy);
279  fscanf(poten_fp, "delta %lf %lf %lf\n", &e.xz, &e.yz, &e.zz);
280 
281  center = origin + e * 0.5
282  * Position(k_nopad[0]-1, k_nopad[1]-1, k_nopad[2]-1);
283 
284  fscanf(poten_fp, "object %*d class gridconnections counts %*lf %*lf %*lf\n");
285  fscanf(poten_fp, "object %*d class array type double rank 0 items %*d data follows\n");
286 
287  // Calculate inverse tensor
288  BigReal det;
289  det = e.xx*(e.yy*e.zz - e.yz*e.zy) - e.xy*(e.yx*e.zz - e.yz*e.zx)
290  + e.xz*(e.yx*e.zy - e.yy*e.zx);
291  inv.xx = (e.yy*e.zz - e.yz*e.zy)/det;
292  inv.xy = -(e.xy*e.zz - e.xz*e.zy)/det;
293  inv.xz = (e.xy*e.yz - e.xz*e.yy)/det;
294  inv.yx = -(e.yx*e.zz - e.yz*e.zx)/det;
295  inv.yy = (e.xx*e.zz - e.xz*e.zx)/det;
296  inv.yz = -(e.xx*e.yz - e.xz*e.yx)/det;
297  inv.zx = (e.yx*e.zy - e.yy*e.zx)/det;
298  inv.zy = -(e.xx*e.zy - e.xy*e.zx)/det;
299  inv.zz = (e.xx*e.yy - e.xy*e.yx)/det;
300 
301  DebugM(4, "origin = " << origin << "\n");
302  DebugM(4, "e = " << e << "\n");
303  DebugM(4, "inv = " << inv << "\n" << endi);
304 }
305 
306 void GridforceFullBaseGrid::readSubgridHierarchy(FILE *poten, int &totalGrids)
307 {
308  DebugM(4, "Beginning of readSubgridHierarchy, generation = " << generation << ", totalGrids = " << totalGrids << "\n" << endi);
309 
310  int elems, generation_in;
311 
313 
314  for (int i = 0; i < numSubgrids; i++) {
315  subgrids[i] = new GridforceFullSubGrid(this);
316  elems = fscanf(poten_fp, "# namdnugrid subgrid %d generation %d min %d %d %d max %d %d %d subgrids count %d\n",
317  &subgrids[i]->subgridIdx, &generation_in,
318  &subgrids[i]->pmin[0], &subgrids[i]->pmin[1], &subgrids[i]->pmin[2],
319  &subgrids[i]->pmax[0], &subgrids[i]->pmax[1], &subgrids[i]->pmax[2],
320  &subgrids[i]->numSubgrids);
321  if (elems < 9) {
322  char msg[256];
323  sprintf(msg, "Problem reading Gridforce potential file! (%d < 9)", elems);
324  NAMD_die(msg);
325  }
326 
327  totalGrids++;
328 
329  if (subgrids[i]->subgridIdx != (totalGrids - 1)) {
330  char msg[256];
331  sprintf(msg, "Problem reading Gridforce potential file! (%d != %d)", subgrids[i]->subgridIdx, totalGrids - 1);
332  NAMD_die(msg);
333  }
334  if (subgrids[i]->generation != generation_in) {
335  char msg[256];
336  sprintf(msg, "Problem reading Gridforce potential file! (%d != %d)", subgrids[i]->generation, generation_in);
337  NAMD_die(msg);
338  }
339 
340 // DebugM(3, "setting maingrid\n");
341 // subgrids[i]->maingrid->subgrids_flat[subgrids[i]->subgridIdx] = subgrids[i];
342 // DebugM(3, "reading subgrid hierarchy\n");
343 
344  subgrids[i]->readSubgridHierarchy(poten, totalGrids);
345  }
346 }
347 
348 
349 /*************************/
350 /* GRIDFORCEFULLMAINGRID */
351 /*************************/
352 
354 {
355  mygridnum = gridnum;
356  generation = 0;
357  subgrids_flat = NULL;
359 }
360 
361 
363 {
364  delete[] subgrids_flat;
365 }
366 
367 
369 {
370  DebugM(4, "Packing maingrid\n" << endi);
371 
372 // msg->put(3*sizeof(float), (char*)pad_p);
373 // msg->put(3*sizeof(float), (char*)pad_n);
374  msg->put(totalGrids);
375  msg->put(mygridnum);
376  msg->put(129*sizeof(char), (char*)filename);
377 
378  DebugM(3, "calling GridforceFullBaseGrid::pack\n" << endi);
379 
381 }
382 
383 
385 {
386  DebugM(4, "Unpacking maingrid\n" << endi);
387 
388 // msg->get(3*sizeof(float), (char*)pad_p);
389 // msg->get(3*sizeof(float), (char*)pad_n);
390  msg->get(totalGrids);
391  msg->get(mygridnum);
392  msg->get(129*sizeof(char), (char*)filename);
393 
395 
396  DebugM(4, "size = " << size << "\n");
397  DebugM(4, "numSubgrids = " << numSubgrids << "\n");
398  DebugM(4, "gapinv = " << gapinv[0] << " " << gapinv[2] << " " << gapinv[2] << " " << "\n");
399  DebugM(4, "generation = " << generation << "\n" << endi);
400  DebugM(4, "filename = " << filename << "\n" << endi);
401 
403 }
404 
405 
407 {
408  DebugM(4, "buildSubgridsFlat() called, totalGrids-1 = " << totalGrids-1 << "\n" << endi);
409  delete[] subgrids_flat;
411  for (int i = 0; i < numSubgrids; i++) {
412  DebugM(3, "adding to subgridsFlat\n" << endi);
414  DebugM(3, "success!\n" << endi);
415  }
416  for (int i = 0; i < totalGrids-1; i++) {
417  DebugM(4, "subgrids_flat[" << i << "]->numSubgrids = " << subgrids_flat[i]->numSubgrids << "\n" << endi);
418  }
419  for (int i = 0; i < numSubgrids; i++) {
420  DebugM(4, "subgrids[" << i << "]->numSubgrids = " << subgrids[i]->numSubgrids << "\n" << endi);
421  }
422 }
423 
424 
425 void GridforceFullMainGrid::initialize(char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams, int brd)
426 {
427  if (brd >= 0) {
428  border = brd;
429  } else {
431  }
432 
433  // FROM init1
434  //FILE *poten = Fopen(potfilename, "r");
435  poten_fp = Fopen(potfilename, "r");
436  if (!poten_fp) {
437  NAMD_die("Problem reading grid force potential file");
438  }
439 
440  // save file name so that grid can be re-read via Tcl
441  strcpy(filename, potfilename);
442 
443  // Read special comment fields and create subgrid objects
444  totalGrids = 1;
445  char line[256];
446  Bool flag = FALSE;
447  numSubgrids = 0;
448  float version;
449  long int poten_offset;
450  do {
451  poten_offset = ftell(poten_fp);
452  fgets(line, 256, poten_fp); // Read comment lines
453  //flag = sscanf(line, "# maingrid subgrids count %d\n", &numSubgrids);
454  flag = sscanf(line, "# namdnugrid version %f\n", &version);
455  } while (line[0] == '#' && !flag);
456 
457  if (flag) {
458  if (version != 1.0) {
459  NAMD_die("Unsupported version of non-uniform grid file format!");
460  }
461  fscanf(poten_fp, "# namdnugrid maingrid subgrids count %d\n", &numSubgrids);
464  } else {
465  fseek(poten_fp, poten_offset, SEEK_SET);
466  }
467 
468  // Read header
469  readHeader(simParams, mgridParams);
470 
471  factor = 1.0;
472  if (mgridParams->gridforceVolts)
473  {
474  factor /= 0.0434; // convert V -> kcal/mol*e
475  }
476  scale = mgridParams->gridforceScale;
477  checksize = mgridParams->gridforceCheckSize;
478 
479  // Allocate storage for potential and read it
480  float *grid_nopad = new float[size_nopad];
481 
482  float tmp2;
483  for (long int count = 0; count < size_nopad; count++) {
484  int err = fscanf(poten_fp, "%f", &tmp2);
485  if (err == EOF || err == 0) {
486  NAMD_die("Grid force potential file incorrectly formatted");
487  }
488  grid_nopad[count] = tmp2 * factor; // temporary, so just store flat
489  }
490  fscanf(poten_fp, "\n");
491 
492  // Shortcuts for accessing 1-D array with four indices
493  dk_nopad[0] = k_nopad[1] * k_nopad[2];
494  dk_nopad[1] = k_nopad[2];
495  dk_nopad[2] = 1;
496 
497  Vector Kvec[3];
498  Kvec[0] = e * Position(k_nopad[0]-1, 0, 0);
499  Kvec[1] = e * Position(0, k_nopad[1]-1, 0);
500  Kvec[2] = e * Position(0, 0, k_nopad[2]-1);
501  Vector Avec[3];
502  Avec[0] = simParams->lattice.a();
503  Avec[1] = simParams->lattice.b();
504  Avec[2] = simParams->lattice.c();
505 
506  // Decide whether we're wrapping
507  for (int i0 = 0; i0 < 3; i0++) {
508  if (mgridParams->gridforceCont[i0])
509  {
510  Bool found = FALSE;
511  for (int i1 = 0; i1 < 3; i1++) {
512  if (cross(Avec[i0].unit(), Kvec[i1].unit()).length() < 1e-4) {
513  found = TRUE;
514  cont[i1] = TRUE;
515  offset[i1] = mgridParams->gridforceVOffset[i0] * factor;
516  // want in grid-point units (normal = 1)
517  gap[i1] = (inv * (Avec[i0] - Kvec[i1])).length();
518  gapinv[i1] = 1.0/gap[i1];
519 
520  if (gap[i1] < 0) {
521  NAMD_die("Gridforce Grid overlap!");
522  }
523 
524  DebugM(4, "cont[" << i1 << "] = " << cont[i1] << "\n");
525  DebugM(4, "gap[" << i1 << "] = " << gap[i1] << "\n");
526  DebugM(4, "gapinv[" << i1 << "] = " << gapinv[i1] << "\n" << endi);
527  }
528  }
529 
530  if (!found) {
531  NAMD_die("No Gridforce unit vector found parallel to requested continuous grid direction!");
532  }
533  } else {
534  // check for grid overlap in non-wrapping dimensions
535  // occurs below
536  }
537  }
538 
539  // Figure out size of true grid (padded on non-periodic sides)
540  Vector delta = 0;
541  for (int i = 0; i < 3; i++) {
542  if (cont[i]) {
543  k[i] = k_nopad[i];
544  } else {
545  k[i] = k_nopad[i] + 2*border;
546  delta[i] -= border;
547  }
548  }
549  DebugM(4, "delta = " << e * delta << " (" << delta << ")\n" << endi);
550  origin += e * delta;
551 
552  // Check for grid overlap
553  if (!fits_lattice(simParams->lattice)) {
554  char errmsg[512];
555  if (checksize) {
556  sprintf(errmsg, "Warning: Periodic cell basis too small for Gridforce grid %d. Set gridforcechecksize off in configuration file to ignore.\n", mygridnum);
557  NAMD_die(errmsg);
558  }
559  }
560 
561  size = k[0] * k[1] * k[2];
562  dk[0] = k[1] * k[2];
563  dk[1] = k[2];
564  dk[2] = 1;
565 
566  DebugM(3, "size = " << size << ", size_nopad = " << size_nopad << "\n" << endi);
567 
568  delete[] grid;
569  grid = new float[size];
570 
571  n_sum[0] = n_sum[1] = n_sum[2] = 0;
572  p_sum[0] = p_sum[1] = p_sum[2] = 0;
573  for (int i0 = 0; i0 < k_nopad[0]; i0++) {
574  for (int i1 = 0; i1 < k_nopad[1]; i1++) {
575  for (int i2 = 0; i2 < k_nopad[2]; i2++) {
576  // Edges are special cases -- take force there to be
577  // zero for smooth transition across potential
578  // boundary
579 
580  long int ind_nopad = i0*dk_nopad[0] + i1*dk_nopad[1] + i2*dk_nopad[2];
581  int j0 = (cont[0]) ? i0 : i0 + border;
582  int j1 = (cont[1]) ? i1 : i1 + border;
583  int j2 = (cont[2]) ? i2 : i2 + border;
584  long int ind = j0*dk[0] + j1*dk[1] + j2*dk[2];
585 
586  if (i0 == 0) n_sum[0] += grid_nopad[ind_nopad];
587  else if (i0 == k_nopad[0]-1) p_sum[0] += grid_nopad[ind_nopad];
588  if (i1 == 0) n_sum[1] += grid_nopad[ind_nopad];
589  else if (i1 == k_nopad[1]-1) p_sum[1] += grid_nopad[ind_nopad];
590  if (i2 == 0) n_sum[2] += grid_nopad[ind_nopad];
591  else if (i2 == k_nopad[2]-1) p_sum[2] += grid_nopad[ind_nopad];
592 
593  //grid[ind] = grid_nopad[ind_nopad];
594  set_grid(j0, j1, j2, grid_nopad[ind_nopad]);
595  }
596  }
597  }
598 
599  const BigReal modThresh = 1.0;
600 
601  BigReal n_avg[3], p_avg[3];
602  int i0;
603  for (int i0 = 0; i0 < 3; i0++) {
604  int i1 = (i0 + 1) % 3;
605  int i2 = (i0 + 2) % 3;
606  n_avg[i0] = n_sum[i0] / (k_nopad[i1] * k_nopad[i2]);
607  p_avg[i0] = p_sum[i0] / (k_nopad[i1] * k_nopad[i2]);
608 
609  if (cont[i0] && fabs(offset[i0] - (p_avg[i0]-n_avg[i0])) > modThresh)
610  {
611  iout << iWARN << "GRID FORCE POTENTIAL DIFFERENCE IN K" << i0
612  << " DIRECTION IS "
613  << offset[i0] - (p_avg[i0]-n_avg[i0])
614  << " KCAL/MOL*E\n" << endi;
615  }
616  }
617 
618  Bool twoPadVals = (cont[0] + cont[1] + cont[2] == 2);
619  double padVal = 0.0;
620  long int weight = 0;
621  if (!twoPadVals) {
622  // Determine pad value (must average)
623  if (!cont[0]) {
624  padVal += p_sum[0] + n_sum[0];
625  weight += 2 * k_nopad[1] * k_nopad[2];
626  }
627  if (!cont[1]) {
628  padVal += p_sum[1] + n_sum[1];
629  weight += 2 * k_nopad[0] * k_nopad[2];
630  }
631  if (!cont[2]) {
632  padVal += p_sum[2] + n_sum[2];
633  weight += 2 * k_nopad[0] * k_nopad[1];
634  }
635  padVal /= weight;
636  }
637 
638  for (int i = 0; i < 3; i++) {
639  pad_n[i] = (cont[i]) ? 0.0 : (twoPadVals) ? n_avg[i] : padVal;
640  pad_p[i] = (cont[i]) ? 0.0 : (twoPadVals) ? p_avg[i] : padVal;
641  DebugM(4, "pad_n[" << i << "] = " << pad_n[i] << "\n");
642  DebugM(4, "pad_p[" << i << "] = " << pad_p[i] << "\n" << endi);
643  }
644 
645  if (cont[0] && cont[1] && cont[2]) {
646  // Nothing to do
647  return;
648  }
649 
650  // Now fill in rest of new grid
651  for (int i0 = 0; i0 < k[0]; i0++) {
652  for (int i1 = 0; i1 < k[1]; i1++) {
653  for (int i2 = 0; i2 < k[2]; i2++) {
654  if ( (cont[0] || (i0 >= border && i0 < k[0]-border))
655  && (cont[1] || (i1 >= border && i1 < k[1]-border))
656  && (cont[2] || i2 == border) )
657  {
658  i2 += k_nopad[2]-1;
659  continue;
660  }
661 
662  long int ind = i0*dk[0] + i1*dk[1] + i2*dk[2];
663 
664  Position pos = e * Position(i0, i1, i2);
665  int var[3] = {i0, i1, i2};
666 
667  for (int dir = 0; dir < 3; dir++) {
668  if (cont[dir])
669  continue;
670 
671  if (var[dir] < border) {
672  //grid[ind] = pad_n[dir];
673  set_grid(i0, i1, i2, pad_n[dir]);
674  } else if (var[dir] >= k[dir]-border) {
675  //grid[ind] = pad_p[dir];
676  set_grid(i0, i1, i2, pad_p[dir]);
677  }
678  }
679 
680 // DebugM(2, "grid[" << ind << "; " << i0 << ", " << i1
681 // << ", " << i2 << "] = " << get_grid(ind)
682 // << "\n" << endi);
683  }
684  }
685  }
686 
687  for (int i0 = 0; i0 < k[0]; i0++) {
688  for (int i1 = 0; i1 < k[1]; i1++) {
689  for (int i2 = 0; i2 < k[2]; i2++) {
690  DebugM(1, "grid[" << i0 << ", " << i1 << ", " << i2 << "] = " << get_grid(i0,i1,i2) << "\n" << endi);
691  }
692  }
693  }
694 
695  // Clean up
696  DebugM(3, "clean up\n" << endi);
697  delete[] grid_nopad;
698 
699  // Call initialize for each subgrid
700  for (int i = 0; i < numSubgrids; i++) {
701  subgrids[i]->poten_fp = poten_fp;
702  subgrids[i]->initialize(simParams, mgridParams);
703  }
704 
705  // close file pointer
706  fclose(poten_fp);
707 }
708 
709 
711 {
712  DebugM(4, "reinitializing grid\n" << endi);
713  initialize(filename, simParams, mgridParams);
714 }
715 
716 
717 long int GridforceFullMainGrid::get_all_gridvals(float **all_gridvals) const
718 {
719  // Creates a flat array of all grid values, including subgrids,
720  // and puts it in the value pointed to by the 'grids'
721  // argument. Returns the resulting array size. Caller is
722  // responsible for destroying the array via 'delete[]'
723 
724  DebugM(4, "get_all_gridvals called\n" << endi);
725 
726  long int sz = 0;
727  sz += size;
728  for (int i = 0; i < totalGrids-1; i++) {
729  sz += subgrids_flat[i]->size;
730  }
731  DebugM(4, "size = " << sz << "\n" << endi);
732 
733  float *grid_vals = new float[sz];
734  long int idx = 0;
735  for (long int i = 0; i < size; i++) {
736  grid_vals[idx++] = grid[i];
737  }
738  for (int j = 0; j < totalGrids-1; j++) {
739  for (long int i = 0; i < subgrids_flat[j]->size; i++) {
740  grid_vals[idx++] = subgrids_flat[j]->grid[i];
741  }
742  }
743  CmiAssert(idx == sz);
744 
745  *all_gridvals = grid_vals;
746 
747  DebugM(4, "get_all_gridvals finished\n" << endi);
748 
749  return sz;
750 }
751 
752 
753 void GridforceFullMainGrid::set_all_gridvals(float *all_gridvals, long int sz)
754 {
755  DebugM(4, "set_all_gridvals called\n" << endi);
756 
757  long int sz_calc = 0;
758  sz_calc += size;
759  for (int i = 0; i < totalGrids-1; i++) {
760  sz_calc += subgrids_flat[i]->size;
761  }
762  CmiAssert(sz == sz_calc);
763 
764  long int idx = 0;
765  for (long int i = 0; i < size; i++) {
766  DebugM(1, "all_gridvals[" << idx << "] = " << all_gridvals[idx] << "\n" << endi);
767  grid[i] = all_gridvals[idx++];
768  }
769  for (int j = 0; j < totalGrids-1; j++) {
770  for (long int i = 0; i < subgrids_flat[j]->size; i++) {
771  DebugM(1, "all_gridvals[" << idx << "] = " << all_gridvals[idx] << "\n" << endi);
772  subgrids_flat[j]->grid[i] = all_gridvals[idx++];
773  }
774  }
775  CmiAssert(idx == sz);
776 
777  DebugM(4, "set_all_gridvals finished\n" << endi);
778 }
779 
780 
781 void GridforceFullMainGrid::compute_b(float *b, int *inds, Vector gapscale) const
782 {
783  for (int i0 = 0; i0 < 8; i0++) {
784  int inds2[3];
785  int zero_derivs = FALSE;
786 
787  float voff = 0.0;
788  int bit = 1; // bit = 2^i1 in the below loop
789  for (int i1 = 0; i1 < 3; i1++) {
790  inds2[i1] = (inds[i1] + ((i0 & bit) ? 1 : 0)) % k[i1];
791 
792  // Deal with voltage offsets
793  if (cont[i1] && inds[i1] == (k[i1]-1) && inds2[i1] == 0) {
794  voff += offset[i1];
795  DebugM(3, "offset[" << i1 << "] = " << offset[i1] << "\n" << endi);
796  }
797 
798  bit <<= 1; // i.e. multiply by 2
799  }
800 
801  DebugM(1, "inds2 = " << inds2[0] << " " << inds2[1] << " " << inds2[2] << "\n" << endi);
802 
803  // NOTE: leaving everything in terms of unit cell coordinates for now,
804  // eventually will multiply by inv tensor when applying the force
805 
806  // First set variables 'dk_{hi,lo}' (glob notation). The 'hi'
807  // ('lo') variable in a given dimension is the number added (subtracted)
808  // to go up (down) one grid point in that dimension; both are normally
809  // just the corresponding 'dk[i]'. However, if we are sitting on a
810  // boundary and we are using a continuous grid, then we want to map the
811  // next point off the grid back around to the other side. e.g. point
812  // (k[0], i1, k) maps to point (0, i1, k), which would be
813  // accomplished by changing 'dk1_hi' to -(k[0]-1)*dk1.
814 
815  int d_hi[3] = {1, 1, 1};
816  int d_lo[3] = {1, 1, 1};
817  float voffs[3];
818  float dscales[3] = {0.5, 0.5, 0.5};
819  for (int i1 = 0; i1 < 3; i1++) {
820  if (inds2[i1] == 0) {
821  if (cont[i1]) {
822  d_lo[i1] = -(k[i1]-1);
823  voffs[i1] = offset[i1];
824  dscales[i1] = 1.0/(1.0 + gap[i1]) * 1.0/gapscale[i1];
825  }
826  else zero_derivs = TRUE;
827  }
828  else if (inds2[i1] == k[i1]-1) {
829  if (cont[i1]) {
830  d_hi[i1] = -(k[i1]-1);
831  voffs[i1] = offset[i1];
832  dscales[i1] = 1.0/(1.0 + gap[i1]) * 1.0/gapscale[i1];
833  }
834  else zero_derivs = TRUE;
835  }
836  else {
837  voffs[i1] = 0.0;
838  }
839  }
840 
841 // DebugM(2, "cont = " << cont[0] << " " << cont[1] << " " << cont[2] << "\n" << endi);
842 // DebugM(2, "zero_derivs = " << zero_derivs << "\n" << endi);
843 // DebugM(2, "d_hi = " << d_hi[0] << " " << d_hi[1] << " " << d_hi[2] << "\n" << endi);
844 // DebugM(2, "d_lo = " << d_lo[0] << " " << d_lo[1] << " " << d_lo[2] << "\n" << endi);
845  DebugM(1, "dscales = " << dscales[0] << " " << dscales[1] << " " << dscales[2] << "\n" << endi);
846  DebugM(1, "voffs = " << voffs[0] << " " << voffs[1] << " " << voffs[2] << "\n" << endi);
847 
848  // V
849  b[i0] = get_grid(inds2[0],inds2[1],inds2[2]) + voff;
850 
851  if (zero_derivs) {
852  DebugM(2, "zero_derivs\n" << endi);
853  b[8+i0] = 0.0;
854  b[16+i0] = 0.0;
855  b[24+i0] = 0.0;
856  b[32+i0] = 0.0;
857  b[40+i0] = 0.0;
858  b[48+i0] = 0.0;
859  b[56+i0] = 0.0;
860  } else {
861  b[8+i0] = dscales[0] * (get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]) + voffs[0]); // dV/dx
862  b[16+i0] = dscales[1] * (get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]) - get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]) + voffs[1]); // dV/dy
863  b[24+i0] = dscales[2] * (get_grid_d(inds2[0],inds2[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0],inds2[1],inds2[2]-d_lo[2]) + voffs[2]); // dV/dz
864  b[32+i0] = dscales[0] * dscales[1]
865  * (get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]) -
866  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2])); // d2V/dxdy
867  b[40+i0] = dscales[0] * dscales[2]
868  * (get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]+d_hi[2]) -
869  get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]-d_lo[2])); // d2V/dxdz
870  b[48+i0] = dscales[1] * dscales[2]
871  * (get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) -
872  get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2])); // d2V/dydz
873 
874  b[56+i0] = dscales[0] * dscales[1] * dscales[2] // d3V/dxdydz
875  * (get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) -
876  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) +
877  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) +
878  get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2]));
879  }
880 
881  DebugM(1, "V = " << b[i0] << "\n");
882 
883  DebugM(1, "dV/dx = " << b[8+i0] << "\n");
884  DebugM(1, "dV/dy = " << b[16+i0] << "\n");
885  DebugM(1, "dV/dz = " << b[24+i0] << "\n");
886 
887  DebugM(1, "d2V/dxdy = " << b[32+i0] << "\n");
888  DebugM(1, "d2V/dxdz = " << b[40+i0] << "\n");
889  DebugM(1, "d2V/dydz = " << b[48+i0] << "\n");
890 
891  DebugM(1, "d3V/dxdydz = " << b[56+i0] << "\n" << endi);
892  }
893 }
894 
895 
896 /************************/
897 /* GRIDFORCEFULLSUBGRID */
898 /************************/
899 
901  parent = parent_in;
905  while (tmp->generation > 0) {
906  tmp = ((GridforceFullSubGrid *)tmp)->parent;
907  }
909  DebugM(4, "generation = " << generation << "\n" << endi);
910 }
911 
912 
914 {
915  int tmp;
916  char line[256];
917  long int poten_offset;
918 
919  // Skip 'attribute's
920  DebugM(3, "Skipping 'attribute' keywords...\n" << endi);
921  char str[256];
922  do {
923  poten_offset = ftell(poten_fp);
924  fscanf(poten_fp, "%s", str);
925  fgets(line, 256, poten_fp);
926  DebugM(4, "Read line " << str << " " << line << endi);
927  } while (strcmp(str, "attribute") == 0);
928  fseek(poten_fp, poten_offset, SEEK_SET);
929 
930  // Skip 'field' object
931  DebugM(3, "Skipping 'field' object\n" << endi);
932  fscanf(poten_fp, "object");
933  int n;
934  n = fscanf(poten_fp, "\"%[^\"]\" class field\n", str);
935  if (n == 0) {
936  n = fscanf(poten_fp, "%d class field\n", &tmp);
937  }
938 
939  if (n == 0) {
940  NAMD_die("Error reading gridforce grid! Could not find field object!\n");
941  }
942 
943  // Skip 'component's
944  DebugM(3, "Skipping 'component' keywords\n" << endi);
945  do {
946  poten_offset = ftell(poten_fp);
947  fscanf(poten_fp, "%s", str);
948  fgets(line, 256, poten_fp);
949  } while (strcmp(str, "component") == 0);
950  fseek(poten_fp, poten_offset, SEEK_SET);
951 
952  // Read header
953  readHeader(simParams, mgridParams);
954 
955  factor = 1.0;
956  if (mgridParams->gridforceVolts)
957  {
958  factor /= 0.0434; // convert V -> kcal/mol*e
959  }
960  scale = mgridParams->gridforceScale;
961  checksize = mgridParams->gridforceCheckSize;
962 
963  for (int i = 0; i < 3; i++) {
964  k[i] = k_nopad[i]; // subgrids aren't padded
965  }
966 
967  // Make sure that each subgrid dimension is an integral
968  // number of spanned supergrid cells. This is to ensure that no
969  // supergrid nodes land in the middle of a subgrid, because in
970  // this case forces could not be matched properly.
971  for (int i = 0; i < 3; i++) {
972  if ((k[i] - 1) % (pmax[i] - pmin[i] + 1) != 0) {
973  iout << (k[i] - 1) << " % " << (pmax[i] - pmin[i] + 1) << " != 0\n" << endi;
974  NAMD_die("Error reading gridforce grid! Subgrid dimensions must be an integral number spanned parent cells!");
975  }
976  }
977 
978  for (int i = 0; i < 3; i++) {
979  if (parent->cont[i]) {
980  cont[i] = (pmin[i] == 0 && pmax[i] == parent->k[i]-2) ? TRUE : FALSE;
981  DebugM(3, "pmin[" << i << "] = " << pmin[i] << " pmax[" << i << "] = " << pmax[i] << " parent->k[" << i << "] = " << parent->k[i] << " cont[" << i << "] = " << cont[i] << "\n" << endi);
982  } else {
983  cont[i] = false;
984  if (parent->generation == 0) {
985  // Add to pmin, pmax since parent got extra gridpoint layer(s) (maybe)
986  int brd = parent->get_border();
987  pmin[i] += brd;
988  pmax[i] += brd;
989  }
990  }
991  }
992 
993  DebugM(4, "pmin = " << pmin[0] << " " << pmin[1] << " " << pmin[2] << "\n");
994  DebugM(4, "pmax = " << pmax[0] << " " << pmax[1] << " " << pmax[2] << "\n" << endi);
995 
996  Vector origin2 = parent->origin + parent->e * Position(pmin[0], pmin[1], pmin[2]);
997  Vector escale, invscale;
998  for (int i = 0; i < 3; i++) {
999  escale[i] = double(pmax[i] - pmin[i] + 1)/(k[i]-1);
1000  invscale[i] = 1.0/escale[i];
1001  if (cont[i]) { pmax[i]++; }
1002  }
1003  Tensor e2 = tensorMult(parent->e, Tensor::diagonal(escale));
1004 
1005  // Check that lattice parameters agree with min and max numbers
1006  // from subgrid hierarchy.
1007  double TOL2 = 1e-4; // Totally arbitrary
1008  if (pow(origin2.x-origin.x, 2) > TOL2 ||
1009  pow(origin2.y-origin.y, 2) > TOL2 ||
1010  pow(origin2.z-origin.z, 2) > TOL2 ||
1011  pow(e2.xx-e.xx, 2) > TOL2 ||
1012  pow(e2.xy-e.xy, 2) > TOL2 ||
1013  pow(e2.xz-e.xz, 2) > TOL2 ||
1014  pow(e2.yx-e.yx, 2) > TOL2 ||
1015  pow(e2.yy-e.yy, 2) > TOL2 ||
1016  pow(e2.yz-e.yz, 2) > TOL2 ||
1017  pow(e2.zx-e.zx, 2) > TOL2 ||
1018  pow(e2.zy-e.zy, 2) > TOL2 ||
1019  pow(e2.zz-e.zz, 2) > TOL2)
1020  {
1021  NAMD_die("Error reading gridforce grid! Subgrid lattice does not match!");
1022  }
1023 
1024  // Overwrite what was read from the header
1025  origin = origin2;
1026  e = e2;
1027 
1028  inv = tensorMult(Tensor::diagonal(invscale), parent->inv);
1029  for (int i = 0; i < 3; i++) {
1030  gap[i] = escale[i] * parent->gap[i];
1031  gapinv[i] = invscale[i] * parent->gapinv[i];
1032  offset[i] = parent->offset[i];
1033  }
1034  center = origin + e * 0.5 * Position(k[0], k[1], k[2]);
1035 
1036  DebugM(4, "origin = " << origin << "\n");
1037  DebugM(4, "e = " << e << "\n");
1038  DebugM(4, "inv = " << inv << "\n");
1039  DebugM(4, "gap = " << gap[0] << " " << gap[2] << " " << gap[2] << " " << "\n");
1040  DebugM(4, "gapinv = " << gapinv[0] << " " << gapinv[2] << " " << gapinv[2] << " " << "\n");
1041  DebugM(4, "numSubgrids = " << numSubgrids << "\n");
1042  DebugM(4, "k = " << k[0] << " " << k[1] << " " << k[2] << "\n");
1043  DebugM(4, "escale = " << escale << "\n");
1044  DebugM(4, "invscale = " << invscale << "\n" << endi);
1045 
1046  /*** Set members ***/
1047  size = k[0] * k[1] * k[2];
1048  dk[0] = k[1] * k[2];
1049  dk[1] = k[2];
1050  dk[2] = 1;
1051 
1052  scale_dV = Tensor::diagonal(escale);
1053  scale_d2V = Tensor::diagonal(Vector(escale.x*escale.y, escale.x*escale.z, escale.y*escale.z));
1054  scale_d3V = escale.x * escale.y * escale.z;
1055 
1056  DebugM(4, "scale_dV = " << scale_dV << "\n");
1057  DebugM(4, "scale_d2V = " << scale_d2V << "\n");
1058  DebugM(4, "scale_d3V = " << scale_d3V << "\n" << endi);
1059 
1060  // Allocate storage for potential and read it
1061  float *grid_tmp = new float[size];
1062 
1063  float tmp2;
1064  DebugM(3, "size_nopad = " << size_nopad << "\n");
1065  for (long int count = 0; count < size_nopad; count++) {
1066 // poten_offset = ftell(poten_fp);
1067 // fscanf(poten_fp, "%s", str);
1068 // fgets(line, 256, poten_fp);
1069 // DebugM(4, "Read line " << str << " " << line << endi);
1070 // fseek(poten_fp, poten_offset, SEEK_SET);
1071 
1072  int err = fscanf(poten_fp, "%f", &tmp2);
1073  if (err == EOF || err == 0) {
1074  NAMD_die("Grid force potential file incorrectly formatted");
1075  }
1076  grid_tmp[count] = tmp2 * factor;
1077  }
1078  fscanf(poten_fp, "\n");
1079 
1080  // Set real grid
1081  DebugM(3, "allocating grid\n" << endi);
1082  delete[] grid;
1083  grid = new float[size];
1084  for (int i0 = 0; i0 < k_nopad[0]; i0++) {
1085  for (int i1 = 0; i1 < k_nopad[1]; i1++) {
1086  for (int i2 = 0; i2 < k_nopad[2]; i2++) {
1087  long int ind = i0*dk[0] + i1*dk[1] + i2*dk[2];
1088  set_grid(i0, i1, i2, grid_tmp[ind]);
1089  }
1090  }
1091  }
1092 
1093  for (int i0 = 0; i0 < k[0]; i0++) {
1094  for (int i1 = 0; i1 < k[1]; i1++) {
1095  for (int i2 = 0; i2 < k[2]; i2++) {
1096  DebugM(1, "grid[" << i0 << ", " << i1 << ", " << i2 << "] = " << get_grid(i0,i1,i2) << "\n" << endi);
1097  }
1098  }
1099  }
1100 
1101  // Clean up
1102  delete[] grid_tmp;
1103 
1104  // Call initialize for each subgrid
1105  for (int i = 0; i < numSubgrids; i++) {
1106  subgrids[i]->initialize(simParams, mgridParams);
1107  }
1108 }
1109 
1110 
1112 {
1113  DebugM(4, "Packing subgrid\n" << endi);
1114 
1115  msg->put(sizeof(Tensor), (char*)&scale_dV);
1116  msg->put(sizeof(Tensor), (char*)&scale_d2V);
1117  msg->put(sizeof(float), (char*)&scale_d3V);
1118 
1119  msg->put(3*sizeof(int), (char*)pmin);
1120  msg->put(3*sizeof(int), (char*)pmax);
1121  msg->put(subgridIdx);
1122 
1123  DebugM(3, "calling GridforceFullBaseGrid::pack\n" << endi);
1124 
1126 }
1127 
1128 
1130 {
1131  DebugM(4, "Unpacking subgrid\n" << endi);
1132 
1133  msg->get(sizeof(Tensor), (char*)&scale_dV);
1134  msg->get(sizeof(Tensor), (char*)&scale_d2V);
1135  msg->get(sizeof(float), (char*)&scale_d3V);
1136 
1137  msg->get(3*sizeof(int), (char*)pmin);
1138  msg->get(3*sizeof(int), (char*)pmax);
1139  msg->get(subgridIdx);
1140 
1142 
1143  DebugM(4, "size = " << size << "\n");
1144  DebugM(4, "numSubgrids = " << numSubgrids << "\n");
1145  DebugM(4, "gapinv = " << gapinv[0] << " " << gapinv[2] << " " << gapinv[2] << " " << "\n");
1146  DebugM(4, "generation = " << generation << "\n" << endi);
1147 }
1148 
1149 
1151 {
1152  DebugM(4, "addToSubgridsFlat() called, subgridIdx = " << subgridIdx << ", maingrid->numSubgrids = " << maingrid->numSubgrids << "\n" << endi);
1153  maingrid->subgrids_flat[subgridIdx-1] = this;
1154  for (int i = 0; i < numSubgrids; i++) {
1156  }
1157 }
1158 
1159 
1160 void GridforceFullSubGrid::compute_b(float *b, int *inds, Vector gapscale) const
1161 {
1162  for (int i0 = 0; i0 < 8; i0++) {
1163  int inds2[3];
1164 
1165  float voff = 0.0;
1166  int bit = 1; // bit = 2^i1 in the below loop
1167  for (int i1 = 0; i1 < 3; i1++) {
1168  inds2[i1] = (inds[i1] + ((i0 & bit) ? 1 : 0)) % k[i1];
1169 
1170  // Deal with voltage offsets
1171  if (cont[i1] && inds[i1] == (k[i1]-1) && inds2[i1] == 0) {
1172  voff += offset[i1];
1173  DebugM(3, "offset[" << i1 << "] = " << offset[i1] << "\n" << endi);
1174  }
1175 
1176  bit <<= 1; // i.e. multiply by 2
1177  }
1178 
1179  DebugM(3, "inds2 = " << inds2[0] << " " << inds2[1] << " " << inds2[2] << "\n" << endi);
1180 
1181  int d_hi[3] = {1, 1, 1};
1182  int d_lo[3] = {1, 1, 1};
1183  float voffs[3];
1184  float dscales[3] = {0.5, 0.5, 0.5};
1185  for (int i1 = 0; i1 < 3; i1++) {
1186  if (inds2[i1] == 0 && cont[i1]) {
1187  d_lo[i1] = -(k[i1]-1);
1188  voffs[i1] = offset[i1];
1189  dscales[i1] = 1.0/(1.0 + gap[i1]) * 1.0/gapscale[i1];
1190  }
1191  else if (inds2[i1] == k[i1]-1 && cont[i1]) {
1192  d_hi[i1] = -(k[i1]-1);
1193  voffs[i1] = offset[i1];
1194  dscales[i1] = 1.0/(1.0 + gap[i1]) * 1.0/gapscale[i1];
1195  }
1196  else {
1197  voffs[i1] = 0.0;
1198  }
1199  }
1200 
1201  bool edge = false;
1202  for (int i1 = 0; i1 < 3; i1++) {
1203  if (!cont[i1] && (inds2[i1] == 0 || inds2[i1] == k[i1]-1)) {
1204  edge = true;
1205  }
1206  }
1207 
1208  if (inds2[2] == 0) {
1209 // DebugM(3, "cont = " << cont[0] << " " << cont[1] << " " << cont[2] << " d_hi = " << d_hi[0] << " " << d_hi[1] << " " << d_hi[2] << " d_lo = " << d_lo[0] << " " << d_lo[1] << " " << d_lo[2] << " dscales = " << dscales[0] << " " << dscales[1] << " " << dscales[2] << "\n" << endi);
1210  DebugM(3, "cont = " << cont[0] << " " << cont[1] << " " << cont[2] << "\n" << endi);
1211  }
1212 
1213  if (edge) {
1214  DebugM(2, "Edge!\n" << endi);
1215 
1216  // Must get derivatives from parent
1217  Position pos = e * Vector(inds2[0], inds2[1], inds2[2]) + origin; // Gridpoint position in realspace
1218  Vector g = parent->inv * (pos - parent->origin); // Gridpoint position in parent's gridspace
1219  Vector dg;
1220  int inds3[3];
1221 
1222  DebugM(2, "g = " << g << "\n" << endi);
1223 
1224  for (int i = 0; i < 3; i++) {
1225  inds3[i] = (int)floor(g[i]);
1226  dg[i] = g[i] - inds3[i];
1227  }
1228 
1229  float x[4], y[4], z[4];
1230  x[0] = 1; y[0] = 1; z[0] = 1;
1231  for (int j = 1; j < 4; j++) {
1232  x[j] = x[j-1] * dg.x;
1233  y[j] = y[j-1] * dg.y;
1234  z[j] = z[j-1] * dg.z;
1235  DebugM(1, "x[" << j << "] = " << x[j] << "\n");
1236  DebugM(1, "y[" << j << "] = " << y[j] << "\n");
1237  DebugM(1, "z[" << j << "] = " << z[j] << "\n" << endi);
1238  }
1239 
1240  // Compute parent matrices
1241  float b_parent[64];
1242  parent->compute_b(b_parent, inds3, gapscale);
1243 
1244  float a_parent[64];
1245  parent->compute_a(a_parent, b_parent);
1246 
1247  // Compute parent derivatives
1248  float V = parent->compute_V(a_parent, x, y, z);
1249  Vector dV = scale_dV * parent->compute_dV(a_parent, x, y, z);
1250  Vector d2V = scale_d2V * parent->compute_d2V(a_parent, x, y, z);
1251  float d3V = scale_d3V * parent->compute_d3V(a_parent, x, y, z);
1252 
1253  b[i0] = V;
1254  b[8+i0] = dV[0];
1255  b[16+i0] = dV[1];
1256  b[24+i0] = dV[2];
1257  b[32+i0] = d2V[0];
1258  b[40+i0] = d2V[1];
1259  b[48+i0] = d2V[2];
1260  b[56+i0] = d3V;
1261  } else {
1262  b[i0] = get_grid(inds2[0],inds2[1],inds2[2]) + voff; // V
1263 
1264  b[8+i0] = dscales[0] * (get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]) + voffs[0]); // dV/dx
1265  b[16+i0] = dscales[1] * (get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]) - get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]) + voffs[1]); // dV/dy
1266  b[24+i0] = dscales[2] * (get_grid_d(inds2[0],inds2[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0],inds2[1],inds2[2]-d_lo[2]) + voffs[2]); // dV/dz
1267  b[32+i0] = dscales[0] * dscales[1]
1268  * (get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]) -
1269  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2])); // d2V/dxdy
1270  b[40+i0] = dscales[0] * dscales[2]
1271  * (get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]+d_hi[2]) -
1272  get_grid_d(inds2[0]+d_hi[0],inds2[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1],inds2[2]-d_lo[2])); // d2V/dxdz
1273  b[48+i0] = dscales[1] * dscales[2]
1274  * (get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) -
1275  get_grid_d(inds2[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2])); // d2V/dydz
1276 
1277  b[56+i0] = dscales[0] * dscales[1] * dscales[2] // d3V/dxdydz
1278  * (get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]+d_hi[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) -
1279  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]+d_hi[2]) +
1280  get_grid_d(inds2[0]+d_hi[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2]) + get_grid_d(inds2[0]-d_lo[0],inds2[1]+d_hi[1],inds2[2]-d_lo[2]) +
1281  get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2]+d_hi[2]) - get_grid_d(inds2[0]-d_lo[0],inds2[1]-d_lo[1],inds2[2]-d_lo[2]));
1282  }
1283 
1284  if (inds2[0] == 1 && inds2[1] == 1 && inds2[2] == 0) {
1285  DebugM(1, "Sub V = " << b[i0] << "\n");
1286 
1287  DebugM(1, "Sub dV/dx = " << b[8+i0] << "\n");
1288  DebugM(1, "Sub dV/dy = " << b[16+i0] << "\n");
1289  DebugM(1, "Sub dV/dz = " << b[24+i0] << "\n");
1290 
1291  DebugM(1, "Sub d2V/dxdy = " << b[32+i0] << "\n");
1292  DebugM(1, "Sub d2V/dxdz = " << b[40+i0] << "\n");
1293  DebugM(1, "Sub d2V/dydz = " << b[48+i0] << "\n");
1294 
1295  DebugM(1, "Sub d3V/dxdydz = " << b[56+i0] << "\n" << endi);
1296  }
1297  }
1298 }
1299 
1300 
1301 /*********************/
1302 /* GRIDFORCELITEGRID */
1303 /*********************/
1304 
1306 {
1307  mygridnum = gridnum;
1308  grid = NULL;
1310 }
1311 
1312 
1314 {
1315  delete[] grid;
1316 }
1317 
1318 
1320 {
1321  // cheat and use GridforceFullMainGrid to read the file
1323  tmp_grid->initialize(potfilename, simParams, mgridParams, 1);
1324 
1325  if (tmp_grid->get_total_grids() != 1) {
1326  NAMD_die("Cannot use gridforcelite option with multi-resolution grid!");
1327  }
1328 
1329  // save file name so that grid can be re-read via Tcl
1330  strcpy(filename, potfilename);
1331 
1332  // copy parameters
1333  k[0] = tmp_grid->get_k0();
1334  k[1] = tmp_grid->get_k1();
1335  k[2] = tmp_grid->get_k2();
1336  k[3] = 4; // for V, dV/dx, dV/dy, dV/dz grids
1337  origin = tmp_grid->get_origin();
1338  center = tmp_grid->get_center();
1339  e = tmp_grid->get_e();
1340  inv = tmp_grid->get_inv();
1341  scale = tmp_grid->get_scale();
1342 
1343  // calculate rest of parameters
1344  size = k[0] * k[1] * k[2] * k[3];
1345  dk[0] = k[1] * k[2] * k[3];
1346  dk[1] = k[2] * k[3];
1347  dk[2] = k[3];
1348  dk[3] = 1;
1349 
1350  // copy the potential grid
1351  delete[] grid;
1352  grid = new float[size];
1353  for (int i0 = 0; i0 < k[0]; i0++) {
1354  for (int i1 = 0; i1 < k[1]; i1++) {
1355  for (int i2 = 0; i2 < k[2]; i2++) {
1356  float V = tmp_grid->get_grid(i0, i1, i2);
1357  set_grid(i0, i1, i2, 0, V);
1358  DebugM(1, "V[" << i0 << "," << i1 << "," << i2 << "] = " << get_grid(i0, i1, i2, 0) << "(" << V << ")\n" << endi);
1359  }
1360  }
1361  }
1362 
1363  delete tmp_grid;
1364 
1366 }
1367 
1368 
1370 {
1371  // calculate derivative grids
1372  // separate loop so all V values have been set already
1373  for (int i0 = 0; i0 < k[0]; i0++) {
1374  for (int i1 = 0; i1 < k[1]; i1++) {
1375  for (int i2 = 0; i2 < k[2]; i2++) {
1376  float dx, dy, dz;
1377  if (i0 == 0 || i0 == k[0]-1 || i1 == 0 || i1 == k[1]-1 || i2 == 0 || i2 == k[2]-1) {
1378  // on edge, set ALL derivatives to zero (make up for lack of padding)
1379  dx = 0;
1380  dy = 0;
1381  dz = 0;
1382  } else {
1383  dx = 0.5 * (get_grid_d(i0+1,i1,i2,0) - get_grid_d(i0-1,i1,i2,0));
1384  dy = 0.5 * (get_grid_d(i0,i1+1,i2,0) - get_grid_d(i0,i1-1,i2,0));
1385  dz = 0.5 * (get_grid_d(i0,i1,i2+1,0) - get_grid_d(i0,i1,i2-1,0));
1386  }
1387  set_grid(i0, i1, i2, 1, dx);
1388  set_grid(i0, i1, i2, 2, dy);
1389  set_grid(i0, i1, i2, 3, dz);
1390  DebugM(1, "dx[" << i0 << "," << i1 << "," << i2 << "] = " << get_grid(i0, i1, i2, 1) << "(" << dx << ")\n" << endi);
1391  DebugM(1, "dy[" << i0 << "," << i1 << "," << i2 << "] = " << get_grid(i0, i1, i2, 2) << "(" << dy << ")\n" << endi);
1392  DebugM(1, "dz[" << i0 << "," << i1 << "," << i2 << "] = " << get_grid(i0, i1, i2, 3) << "(" << dz << ")\n" << endi);
1393  }
1394  }
1395  }
1396 }
1397 
1398 
1400 {
1401  initialize(filename, simParams, mgridParams);
1402 }
1403 
1404 
1406 {
1407  msg->put(4*sizeof(int), (char*)k);
1408  msg->put(size);
1409  msg->put(4*sizeof(long int), (char*)dk);
1410 
1411  msg->put(sizeof(Vector), (char*)&origin);
1412  msg->put(sizeof(Vector), (char*)&center);
1413  msg->put(sizeof(Tensor), (char*)&e);
1414  msg->put(sizeof(Tensor), (char*)&inv);
1415  msg->put(sizeof(Vector), (char*)&scale);
1416  msg->put(sizeof(Bool), (char*)&checksize);
1417 
1418  msg->put(129*sizeof(char), (char*)filename);
1419 
1420  msg->put(size*sizeof(float), (char*)grid);
1421 }
1422 
1423 
1425 {
1426  delete[] grid;
1427  grid = NULL;
1428 
1429  msg->get(4*sizeof(int), (char*)k);
1430  msg->get(size);
1431  msg->get(4*sizeof(long int), (char*)dk);
1432 
1433  msg->get(sizeof(Vector), (char*)&origin);
1434  msg->get(sizeof(Vector), (char*)&center);
1435  msg->get(sizeof(Tensor), (char*)&e);
1436  msg->get(sizeof(Tensor), (char*)&inv);
1437  msg->get(sizeof(Vector), (char*)&scale);
1438  msg->get(sizeof(Bool), (char*)&checksize);
1439 
1440  msg->get(129*sizeof(char), (char*)filename);
1441 
1442  if (size) {
1443  grid = new float[size];
1444  msg->get(size*sizeof(float), (char*)grid);
1445  }
1446 }
1447 
1448 
1449 long int GridforceLiteGrid::get_all_gridvals(float** all_gridvals) const
1450 {
1451  // Creates a flat array of all grid values and puts it in the
1452  // value pointed to by the 'all_gridvals' argument. Returns the
1453  // resulting array size. Caller is responsible for destroying the
1454  // array via 'delete[]'
1455 
1456  DebugM(4, "GridforceLiteGrid::get_all_gridvals called\n" << endi);
1457 
1458  long int sz = size;
1459  DebugM(4, "size = " << sz << "\n" << endi);
1460 
1461  float *grid_vals = new float[sz];
1462  long int idx = 0;
1463  for (long int i = 0; i < size; i++) {
1464  grid_vals[idx++] = grid[i];
1465  }
1466  CmiAssert(idx == sz);
1467 
1468  *all_gridvals = grid_vals;
1469 
1470  DebugM(4, "GridforceLiteGrid::get_all_gridvals finished\n" << endi);
1471 
1472  return sz;
1473 }
1474 
1475 
1476 void GridforceLiteGrid::set_all_gridvals(float* all_gridvals, long int sz)
1477 {
1478  DebugM(4, "GridforceLiteGrid::set_all_gridvals called\n" << endi);
1479 
1480  long int sz_calc = size;
1481  CmiAssert(sz == sz_calc);
1482 
1483  long int idx = 0;
1484  for (long int i = 0; i < size; i++) {
1485  grid[i] = all_gridvals[idx++];
1486  }
1487  CmiAssert(idx == sz);
1488 
1489  //compute_derivative_grids(); // not needed if we're sending all 4 grids
1490 
1491  DebugM(4, "GridforceLiteGrid::set_all_gridvals finished\n" << endi);
1492 }
void reinitialize(SimParameters *simParams, MGridforceParams *mgridParams)
BigReal zy
Definition: Tensor.h:19
void pack(MOStream *msg) const
int get_total_grids(void) const
BigReal xz
Definition: Tensor.h:17
friend class GridforceFullSubGrid
Definition: GridForceGrid.h:83
virtual ~GridforceGrid()
Definition: GridForceGrid.C:48
int get_k2(void) const
float compute_V(float *a, float *x, float *y, float *z) const
static GridforceGrid * new_grid(int gridnum, char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams)
Definition: GridForceGrid.C:34
Definition: Vector.h:72
virtual void compute_b(float *b, int *inds, Vector gapscale) const =0
void initialize(char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams, int border)
long int get_all_gridvals(float **all_gridvals) const
static void pack_grid(GridforceGrid *grid, MOStream *msg)
Definition: GridForceGrid.C:50
Vector compute_d2V(float *a, float *x, float *y, float *z) const
virtual Position get_origin(void) const =0
#define DebugM(x, y)
Definition: Debug.h:75
std::ostream & endi(std::ostream &s)
Definition: InfoStream.C:54
BigReal z
Definition: Vector.h:74
#define FALSE
Definition: common.h:127
int get_k0(void) const
BigReal yz
Definition: Tensor.h:18
double get_grid_d(int i0, int i1, int i2, int i3) const
float compute_d3V(float *a, float *x, float *y, float *z) const
std::ostream & iWARN(std::ostream &s)
Definition: InfoStream.C:82
MIStream * get(char &data)
Definition: MStream.h:29
virtual ~GridforceLiteGrid()
#define iout
Definition: InfoStream.h:51
virtual int get_k2(void) const =0
char filename[NAMD_FILENAME_BUFFER_SIZE]
GridforceFullMainGrid * maingrid
void readHeader(SimParameters *simParams, MGridforceParams *mgridParams)
void unpack(MIStream *msg)
void compute_b(float *b, int *inds, Vector gapscale) const
void compute_derivative_grids(void)
void set_grid(int i0, int i1, int i2, int i3, float V)
GridforceFullSubGrid ** subgrids
GridforceFullBaseGrid * parent
bool fits_lattice(const Lattice &lattice)
Definition: GridForceGrid.C:84
virtual void pack(MOStream *msg) const
Tensor get_e(void) const
void set_grid(int i0, int i1, int i2, float V)
virtual void pack(MOStream *msg) const =0
void NAMD_bug(const char *err_msg)
Definition: common.C:195
BigReal yx
Definition: Tensor.h:18
void compute_b(float *b, int *inds, Vector gapscale) const
Position get_center(void) const
GridforceFullSubGrid ** subgrids_flat
GridforceFullMainGrid(int gridnum)
virtual Position get_center(void) const =0
void reinitialize(SimParameters *simParams, MGridforceParams *mgridParams)
static NAMD_HOST_DEVICE Tensor diagonal(const Vector &v1)
Definition: Tensor.h:37
virtual ~GridforceFullBaseGrid()
int Bool
Definition: common.h:142
void set_all_gridvals(float *all_gridvals, long int sz)
FILE * Fopen(const char *filename, const char *mode)
Definition: common.C:341
long int get_all_gridvals(float **all_gridvals) const
float get_grid(int i0, int i1, int i2) const
Vector Position
Definition: NamdTypes.h:24
GridforceFullSubGrid(GridforceFullBaseGrid *parent_in)
virtual void initialize(char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams)=0
BigReal x
Definition: Vector.h:74
void readSubgridHierarchy(FILE *poten, int &totalGrids)
GridforceLiteGrid(int gridnum)
void addToSubgridsFlat(void)
void NAMD_die(const char *err_msg)
Definition: common.C:147
virtual int get_border(void) const =0
virtual void unpack(MIStream *msg)
void set_all_gridvals(float *all_gridvals, long int sz)
Vector get_scale(void) const
BigReal xx
Definition: Tensor.h:17
BigReal zz
Definition: Tensor.h:19
GridforceGridType type
Definition: GridForceGrid.h:73
char filename[NAMD_FILENAME_BUFFER_SIZE]
Tensor get_inv(void) const
#define simParams
Definition: Output.C:129
void initialize(char *potfilename, SimParameters *simParams, MGridforceParams *mgridParams)
static GridforceGrid * unpack_grid(int gridnum, MIStream *msg)
Definition: GridForceGrid.C:60
void buildSubgridsFlat(void)
Definition: Tensor.h:15
BigReal xy
Definition: Tensor.h:17
Tensor tensorMult(const Tensor &t1, const Tensor &t2)
virtual int get_k0(void) const =0
BigReal y
Definition: Vector.h:74
Position get_corner(int idx)
virtual int get_k1(void) const =0
Position wrap_position(const Position &pos, const Lattice &lattice)
virtual void unpack(MIStream *msg)=0
virtual Tensor get_e(void) const =0
Vector compute_dV(float *a, float *x, float *y, float *z) const
BigReal yy
Definition: Tensor.h:18
void pack(MOStream *msg) const
virtual ~GridforceFullMainGrid()
MOStream * put(char data)
Definition: MStream.h:112
GridforceGridType get_grid_type(void)
Definition: GridForceGrid.h:64
static const int default_border
float get_grid(int i0, int i1, int i2, int i3) const
void pack(MOStream *msg) const
Position get_origin(void) const
double get_grid_d(int i0, int i1, int i2) const
void initialize(SimParameters *simParams, MGridforceParams *mgridParams)
int get_k1(void) const
BigReal zx
Definition: Tensor.h:19
#define TRUE
Definition: common.h:128
double BigReal
Definition: common.h:123
void unpack(MIStream *msg)
void unpack(MIStream *msg)
void compute_a(float *a, float *b) const