1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642
|
/* ----------------------------------------------------------------------
This is the
██╗ ██╗ ██████╗ ██████╗ ██████╗ ██╗ ██╗████████╗███████╗
██║ ██║██╔════╝ ██╔════╝ ██╔════╝ ██║ ██║╚══██╔══╝██╔════╝
██║ ██║██║ ███╗██║ ███╗██║ ███╗███████║ ██║ ███████╗
██║ ██║██║ ██║██║ ██║██║ ██║██╔══██║ ██║ ╚════██║
███████╗██║╚██████╔╝╚██████╔╝╚██████╔╝██║ ██║ ██║ ███████║
╚══════╝╚═╝ ╚═════╝ ╚═════╝ ╚═════╝ ╚═╝ ╚═╝ ╚═╝ ╚══════╝®
DEM simulation engine, released by
DCS Computing Gmbh, Linz, Austria
http://www.dcs-computing.com, office@dcs-computing.com
LIGGGHTS® is part of CFDEM®project:
http://www.liggghts.com | http://www.cfdem.com
Core developer and main author:
Christoph Kloss, christoph.kloss@dcs-computing.com
LIGGGHTS® is open-source, distributed under the terms of the GNU Public
License, version 2 or later. It is distributed in the hope that it will
be useful, but WITHOUT ANY WARRANTY; without even the implied warranty
of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. You should have
received a copy of the GNU General Public License along with LIGGGHTS®.
If not, see http://www.gnu.org/licenses . See also top-level README
and LICENSE files.
LIGGGHTS® and CFDEM® are registered trade marks of DCS Computing GmbH,
the producer of the LIGGGHTS® software and the CFDEM®coupling software
See http://www.cfdem.com/terms-trademark-policy for details.
-------------------------------------------------------------------------
Contributing author and copyright for this file:
This file is from LAMMPS
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
------------------------------------------------------------------------- */
/* ----------------------------------------------------------------------
Contributing authors:
Benoit Leblanc, Dave Rigby, Paul Saxe (Materials Design)
Reese Jones (Sandia)
------------------------------------------------------------------------- */
#include <stdlib.h>
#include <string.h>
#include "unistd.h"
#include "fix_ave_correlate.h"
#include "update.h"
#include "modify.h"
#include "compute.h"
#include "input.h"
#include "variable.h"
#include "memory.h"
#include "error.h"
#include "force.h"
using namespace LAMMPS_NS;
using namespace FixConst;
enum{COMPUTE,FIX,VARIABLE};
enum{ONE,RUNNING};
enum{AUTO,UPPER,LOWER,AUTOUPPER,AUTOLOWER,FULL};
/* ---------------------------------------------------------------------- */
FixAveCorrelate::FixAveCorrelate(LAMMPS * lmp, int narg, char **arg):
Fix (lmp, narg, arg)
{
if (narg < 7) error->all(FLERR,"Illegal fix ave/correlate command");
MPI_Comm_rank(world,&me);
nevery = force->inumeric(FLERR,arg[3]);
nrepeat = force->inumeric(FLERR,arg[4]);
nfreq = force->inumeric(FLERR,arg[5]);
global_freq = nfreq;
// parse values until one isn't recognized
which = new int[narg-6];
argindex = new int[narg-6];
ids = new char*[narg-6];
value2index = new int[narg-6];
nvalues = 0;
int iarg = 6;
while (iarg < narg) {
if (strncmp(arg[iarg],"c_",2) == 0 ||
strncmp(arg[iarg],"f_",2) == 0 ||
strncmp(arg[iarg],"v_",2) == 0) {
if (arg[iarg][0] == 'c') which[nvalues] = COMPUTE;
else if (arg[iarg][0] == 'f') which[nvalues] = FIX;
else if (arg[iarg][0] == 'v') which[nvalues] = VARIABLE;
int n = strlen(arg[iarg]);
char *suffix = new char[n];
strcpy(suffix,&arg[iarg][2]);
char *ptr = strchr(suffix,'[');
if (ptr) {
if (suffix[strlen(suffix)-1] != ']')
error->all(FLERR,"Illegal fix ave/correlate command");
argindex[nvalues] = atoi(ptr+1);
*ptr = '\0';
} else argindex[nvalues] = 0;
n = strlen(suffix) + 1;
ids[nvalues] = new char[n];
strcpy(ids[nvalues],suffix);
delete [] suffix;
nvalues++;
iarg++;
} else break;
}
// optional args
type = AUTO;
ave = ONE;
startstep = 0;
prefactor = 1.0;
fp = NULL;
overwrite = 0;
char *title1 = NULL;
char *title2 = NULL;
char *title3 = NULL;
while (iarg < narg) {
if (strcmp(arg[iarg],"type") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
if (strcmp(arg[iarg+1],"auto") == 0) type = AUTO;
else if (strcmp(arg[iarg+1],"upper") == 0) type = UPPER;
else if (strcmp(arg[iarg+1],"lower") == 0) type = LOWER;
else if (strcmp(arg[iarg+1],"auto/upper") == 0) type = AUTOUPPER;
else if (strcmp(arg[iarg+1],"auto/lower") == 0) type = AUTOLOWER;
else if (strcmp(arg[iarg+1],"full") == 0) type = FULL;
else error->all(FLERR,"Illegal fix ave/correlate command");
iarg += 2;
} else if (strcmp(arg[iarg],"ave") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
if (strcmp(arg[iarg+1],"one") == 0) ave = ONE;
else if (strcmp(arg[iarg+1],"running") == 0) ave = RUNNING;
else error->all(FLERR,"Illegal fix ave/correlate command");
iarg += 2;
} else if (strcmp(arg[iarg],"start") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
startstep = force->inumeric(FLERR,arg[iarg+1]);
iarg += 2;
} else if (strcmp(arg[iarg],"prefactor") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
prefactor = force->numeric(FLERR,arg[iarg+1]);
iarg += 2;
} else if (strcmp(arg[iarg],"file") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
if (me == 0) {
fp = fopen(arg[iarg+1],"w");
if (fp == NULL) {
char str[512];
sprintf(str,"Cannot open fix ave/correlate file %s",arg[iarg+1]);
error->one(FLERR,str);
}
}
iarg += 2;
} else if (strcmp(arg[iarg],"overwrite") == 0) {
overwrite = 1;
iarg += 1;
} else if (strcmp(arg[iarg],"title1") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
delete [] title1;
int n = strlen(arg[iarg+1]) + 1;
title1 = new char[n];
strcpy(title1,arg[iarg+1]);
iarg += 2;
} else if (strcmp(arg[iarg],"title2") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
delete [] title2;
int n = strlen(arg[iarg+1]) + 1;
title2 = new char[n];
strcpy(title2,arg[iarg+1]);
iarg += 2;
} else if (strcmp(arg[iarg],"title3") == 0) {
if (iarg+2 > narg) error->all(FLERR,"Illegal fix ave/correlate command");
delete [] title3;
int n = strlen(arg[iarg+1]) + 1;
title3 = new char[n];
strcpy(title3,arg[iarg+1]);
iarg += 2;
} else error->all(FLERR,"Illegal fix ave/correlate command");
}
// setup and error check
// for fix inputs, check that fix frequency is acceptable
if (nevery <= 0 || nrepeat <= 0 || nfreq <= 0)
error->all(FLERR,"Illegal fix ave/correlate command");
if (nfreq % nevery)
error->all(FLERR,"Illegal fix ave/correlate command");
if (ave == ONE && nfreq < (nrepeat-1)*nevery)
error->all(FLERR,"Illegal fix ave/correlate command");
if (ave != RUNNING && overwrite)
error->all(FLERR,"Illegal fix ave/correlate command");
for (int i = 0; i < nvalues; i++) {
if (which[i] == COMPUTE) {
int icompute = modify->find_compute(ids[i]);
if (icompute < 0)
error->all(FLERR,"Compute ID for fix ave/correlate does not exist");
if (argindex[i] == 0 && modify->compute[icompute]->scalar_flag == 0)
error->all(FLERR,
"Fix ave/correlate compute does not calculate a scalar");
if (argindex[i] && modify->compute[icompute]->vector_flag == 0)
error->all(FLERR,
"Fix ave/correlate compute does not calculate a vector");
if (argindex[i] && argindex[i] > modify->compute[icompute]->size_vector)
error->all(FLERR,"Fix ave/correlate compute vector "
"is accessed out-of-range");
} else if (which[i] == FIX) {
int ifix = modify->find_fix(ids[i]);
if (ifix < 0)
error->all(FLERR,"Fix ID for fix ave/correlate does not exist");
if (argindex[i] == 0 && modify->fix[ifix]->scalar_flag == 0)
error->all(FLERR,"Fix ave/correlate fix does not calculate a scalar");
if (argindex[i] && modify->fix[ifix]->vector_flag == 0)
error->all(FLERR,"Fix ave/correlate fix does not calculate a vector");
if (argindex[i] && argindex[i] > modify->fix[ifix]->size_vector)
error->all(FLERR,
"Fix ave/correlate fix vector is accessed out-of-range");
if (nevery % modify->fix[ifix]->global_freq)
error->all(FLERR,"Fix for fix ave/correlate "
"not computed at compatible time");
} else if (which[i] == VARIABLE) {
int ivariable = input->variable->find(ids[i]);
if (ivariable < 0)
error->all(FLERR,"Variable name for fix ave/correlate does not exist");
if (input->variable->equalstyle(ivariable) == 0)
error->all(FLERR,
"Fix ave/correlate variable is not equal-style variable");
}
}
// npair = # of correlation pairs to calculate
if (type == AUTO) npair = nvalues;
if (type == UPPER || type == LOWER) npair = nvalues*(nvalues-1)/2;
if (type == AUTOUPPER || type == AUTOLOWER) npair = nvalues*(nvalues+1)/2;
if (type == FULL) npair = nvalues*nvalues;
// print file comment lines
if (fp && me == 0) {
if (title1) fprintf(fp,"%s\n",title1);
else fprintf(fp,"# Time-correlated data for fix %s\n",id);
if (title2) fprintf(fp,"%s\n",title2);
else fprintf(fp,"# Timestep Number-of-time-windows\n");
if (title3) fprintf(fp,"%s\n",title3);
else {
fprintf(fp,"# Index TimeDelta Ncount");
if (type == AUTO)
for (int i = 0; i < nvalues; i++)
fprintf(fp," %s*%s",arg[6+i],arg[6+i]);
else if (type == UPPER)
for (int i = 0; i < nvalues; i++)
for (int j = i+1; j < nvalues; j++)
fprintf(fp," %s*%s",arg[6+i],arg[6+j]);
else if (type == LOWER)
for (int i = 0; i < nvalues; i++)
for (int j = 0; j < i-1; j++)
fprintf(fp," %s*%s",arg[6+i],arg[6+j]);
else if (type == AUTOUPPER)
for (int i = 0; i < nvalues; i++)
for (int j = i; j < nvalues; j++)
fprintf(fp," %s*%s",arg[6+i],arg[6+j]);
else if (type == AUTOLOWER)
for (int i = 0; i < nvalues; i++)
for (int j = 0; j < i; j++)
fprintf(fp," %s*%s",arg[6+i],arg[6+j]);
else if (type == FULL)
for (int i = 0; i < nvalues; i++)
for (int j = 0; j < nvalues; j++)
fprintf(fp," %s*%s",arg[6+i],arg[6+j]);
fprintf(fp,"\n");
}
filepos = ftell(fp);
}
delete [] title1;
delete [] title2;
delete [] title3;
// allocate and initialize memory for averaging
// set count and corr to zero since they accumulate
// also set save versions to zero in case accessed via compute_array()
memory->create(values,nrepeat,nvalues,"ave/correlate:values");
memory->create(count,nrepeat,"ave/correlate:count");
memory->create(save_count,nrepeat,"ave/correlate:save_count");
memory->create(corr,nrepeat,npair,"ave/correlate:corr");
memory->create(save_corr,nrepeat,npair,"ave/correlate:save_corr");
int i,j;
for (i = 0; i < nrepeat; i++) {
save_count[i] = count[i] = 0;
for (j = 0; j < npair; j++)
save_corr[i][j] = corr[i][j] = 0.0;
}
// this fix produces a global array
array_flag = 1;
size_array_rows = nrepeat;
size_array_cols = npair+2;
extarray = 0;
// nvalid = next step on which end_of_step does something
// add nvalid to all computes that store invocation times
// since don't know a priori which are invoked by this fix
// once in end_of_step() can set timestep for ones actually invoked
lastindex = -1;
firstindex = 0;
nsample = 0;
nvalid = nextvalid();
modify->addstep_compute_all(nvalid);
}
/* ---------------------------------------------------------------------- */
FixAveCorrelate::~FixAveCorrelate()
{
delete [] which;
delete [] argindex;
delete [] value2index;
for (int i = 0; i < nvalues; i++) delete [] ids[i];
delete [] ids;
memory->destroy(values);
memory->destroy(count);
memory->destroy(save_count);
memory->destroy(corr);
memory->destroy(save_corr);
if (fp && me == 0) fclose(fp);
}
/* ---------------------------------------------------------------------- */
int FixAveCorrelate::setmask()
{
int mask = 0;
mask |= END_OF_STEP;
return mask;
}
/* ---------------------------------------------------------------------- */
void FixAveCorrelate::init()
{
// set current indices for all computes,fixes,variables
for (int i = 0; i < nvalues; i++) {
if (which[i] == COMPUTE) {
int icompute = modify->find_compute(ids[i]);
if (icompute < 0)
error->all(FLERR,"Compute ID for fix ave/correlate does not exist");
value2index[i] = icompute;
} else if (which[i] == FIX) {
int ifix = modify->find_fix(ids[i]);
if (ifix < 0)
error->all(FLERR,"Fix ID for fix ave/correlate does not exist");
value2index[i] = ifix;
} else if (which[i] == VARIABLE) {
int ivariable = input->variable->find(ids[i]);
if (ivariable < 0)
error->all(FLERR,"Variable name for fix ave/correlate does not exist");
value2index[i] = ivariable;
}
}
// need to reset nvalid if nvalid < ntimestep b/c minimize was performed
if (nvalid < update->ntimestep) {
lastindex = -1;
firstindex = 0;
nsample = 0;
nvalid = nextvalid();
modify->addstep_compute_all(nvalid);
}
}
/* ----------------------------------------------------------------------
only does something if nvalid = current timestep
------------------------------------------------------------------------- */
void FixAveCorrelate::setup(int vflag)
{
end_of_step();
}
/* ---------------------------------------------------------------------- */
void FixAveCorrelate::end_of_step()
{
int i,j,m;
double scalar = 0.0;
// skip if not step which requires doing something
bigint ntimestep = update->ntimestep;
if (ntimestep != nvalid) return;
// accumulate results of computes,fixes,variables to origin
// compute/fix/variable may invoke computes so wrap with clear/add
modify->clearstep_compute();
// lastindex = index in values ring of latest time sample
lastindex++;
if (lastindex == nrepeat) lastindex = 0;
for (i = 0; i < nvalues; i++) {
m = value2index[i];
// invoke compute if not previously invoked
if (which[i] == COMPUTE) {
Compute *compute = modify->compute[m];
if (argindex[i] == 0) {
if (!(compute->invoked_flag & INVOKED_SCALAR)) {
compute->compute_scalar();
compute->invoked_flag |= INVOKED_SCALAR;
}
scalar = compute->scalar;
} else {
if (!(compute->invoked_flag & INVOKED_VECTOR)) {
compute->compute_vector();
compute->invoked_flag |= INVOKED_VECTOR;
}
scalar = compute->vector[argindex[i]-1];
}
// access fix fields, guaranteed to be ready
} else if (which[i] == FIX) {
if (argindex[i] == 0)
scalar = modify->fix[m]->compute_scalar();
else
scalar = modify->fix[m]->compute_vector(argindex[i]-1);
// evaluate equal-style variable
} else if (which[i] == VARIABLE)
scalar = input->variable->compute_equal(m);
values[lastindex][i] = scalar;
}
// fistindex = index in values ring of earliest time sample
// nsample = number of time samples in values ring
if (nsample < nrepeat) nsample++;
else {
firstindex++;
if (firstindex == nrepeat) firstindex = 0;
}
nvalid += nevery;
modify->addstep_compute(nvalid);
// calculate all Cij() enabled by latest values
accumulate();
if (ntimestep % nfreq) return;
// save results in save_count and save_corr
for (i = 0; i < nrepeat; i++) {
save_count[i] = count[i];
if (count[i])
for (j = 0; j < npair; j++)
save_corr[i][j] = prefactor*corr[i][j]/count[i];
else
for (j = 0; j < npair; j++)
save_corr[i][j] = 0.0;
}
// output result to file
if (fp && me == 0) {
if (overwrite) fseek(fp,filepos,SEEK_SET);
fprintf(fp,BIGINT_FORMAT " %d\n",ntimestep,nrepeat);
for (i = 0; i < nrepeat; i++) {
fprintf(fp,"%d %d %d",i+1,i*nevery,count[i]);
if (count[i])
for (j = 0; j < npair; j++)
fprintf(fp," %g",prefactor*corr[i][j]/count[i]);
else
for (j = 0; j < npair; j++)
fprintf(fp," 0.0");
fprintf(fp,"\n");
}
fflush(fp);
if (overwrite) {
long fileend = ftell(fp);
ftruncate(fileno(fp),fileend);
}
}
// zero accumulation if requested
// recalculate Cij(0)
if (ave == ONE) {
for (i = 0; i < nrepeat; i++) {
count[i] = 0;
for (j = 0; j < npair; j++)
corr[i][j] = 0.0;
}
nsample = 1;
accumulate();
}
}
/* ----------------------------------------------------------------------
accumulate correlation data using more recently added values
------------------------------------------------------------------------- */
void FixAveCorrelate::accumulate()
{
int i,j,k,m,n,ipair;
for (k = 0; k < nsample; k++) count[k]++;
if (type == AUTO) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++) {
corr[k][ipair++] += values[m][i]*values[n][i];
}
m--;
if (m < 0) m = nrepeat-1;
}
} else if (type == UPPER) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++)
for (j = i+1; j < nvalues; j++)
corr[k][ipair++] += values[m][i]*values[n][j];
m--;
if (m < 0) m = nrepeat-1;
}
} else if (type == LOWER) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++)
for (j = 0; j < i; j++)
corr[k][ipair++] += values[m][i]*values[n][j];
m--;
if (m < 0) m = nrepeat-1;
}
} else if (type == AUTOUPPER) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++)
for (j = i; j < nvalues; j++)
corr[k][ipair++] += values[m][i]*values[n][j];
m--;
if (m < 0) m = nrepeat-1;
}
} else if (type == AUTOLOWER) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++)
for (j = 0; j <= i; j++)
corr[k][ipair++] += values[m][i]*values[n][j];
m--;
if (m < 0) m = nrepeat-1;
}
} else if (type == FULL) {
m = n = lastindex;
for (k = 0; k < nsample; k++) {
ipair = 0;
for (i = 0; i < nvalues; i++)
for (j = 0; j < nvalues; j++)
corr[k][ipair++] += values[m][i]*values[n][j];
m--;
if (m < 0) m = nrepeat-1;
}
}
}
/* ----------------------------------------------------------------------
return I,J array value
------------------------------------------------------------------------- */
double FixAveCorrelate::compute_array(int i, int j)
{
if (j == 0) return 1.0*i*nevery;
else if (j == 1) return 1.0*save_count[i];
else if (save_count[i]) return save_corr[i][j-2];
return 0.0;
}
/* ----------------------------------------------------------------------
nvalid = next step on which end_of_step does something
this step if multiple of nevery, else next multiple
startstep is lower bound
------------------------------------------------------------------------- */
bigint FixAveCorrelate::nextvalid()
{
bigint nvalid = update->ntimestep;
if (startstep > nvalid) nvalid = startstep;
if (nvalid % nevery) nvalid = (nvalid/nevery)*nevery + nevery;
return nvalid;
}
/* ---------------------------------------------------------------------- */
void FixAveCorrelate::reset_timestep(bigint ntimestep)
{
if (ntimestep > nvalid) error->all(FLERR,"Fix ave/correlate missed timestep");
}
|