avoid overflows when computing memory offsets and allocating memory
This commit is contained in:
@ -597,7 +597,7 @@ struct remap_plan_3d *remap_3d_create_plan(
|
|||||||
if (memory == 1) {
|
if (memory == 1) {
|
||||||
if (nrecv > 0) {
|
if (nrecv > 0) {
|
||||||
plan->scratch =
|
plan->scratch =
|
||||||
(FFT_SCALAR *) malloc(nqty*out.isize*out.jsize*out.ksize *
|
(FFT_SCALAR *) malloc((size_t)nqty*out.isize*out.jsize*out.ksize *
|
||||||
sizeof(FFT_SCALAR));
|
sizeof(FFT_SCALAR));
|
||||||
if (plan->scratch == nullptr) return nullptr;
|
if (plan->scratch == nullptr) return nullptr;
|
||||||
}
|
}
|
||||||
|
|||||||
@ -112,7 +112,7 @@ void PairEAMOpt::eval()
|
|||||||
int ntypes2 = ntypes*ntypes;
|
int ntypes2 = ntypes*ntypes;
|
||||||
|
|
||||||
fast_alpha_t* _noalias fast_alpha =
|
fast_alpha_t* _noalias fast_alpha =
|
||||||
(fast_alpha_t*) malloc(ntypes2*(nr+1)*sizeof(fast_alpha_t));
|
(fast_alpha_t*) malloc((size_t)ntypes2*(nr+1)*sizeof(fast_alpha_t));
|
||||||
for (i = 0; i < ntypes; i++) for (j = 0; j < ntypes; j++) {
|
for (i = 0; i < ntypes; i++) for (j = 0; j < ntypes; j++) {
|
||||||
fast_alpha_t* _noalias tab = &fast_alpha[i*ntypes*nr+j*nr];
|
fast_alpha_t* _noalias tab = &fast_alpha[i*ntypes*nr+j*nr];
|
||||||
if (type2rhor[i+1][j+1] >= 0) {
|
if (type2rhor[i+1][j+1] >= 0) {
|
||||||
@ -135,7 +135,7 @@ void PairEAMOpt::eval()
|
|||||||
fast_alpha_t* _noalias tabeight = fast_alpha;
|
fast_alpha_t* _noalias tabeight = fast_alpha;
|
||||||
|
|
||||||
fast_gamma_t* _noalias fast_gamma =
|
fast_gamma_t* _noalias fast_gamma =
|
||||||
(fast_gamma_t*) malloc(ntypes2*(nr+1)*sizeof(fast_gamma_t));
|
(fast_gamma_t*) malloc((size_t)ntypes2*(nr+1)*sizeof(fast_gamma_t));
|
||||||
for (i = 0; i < ntypes; i++) for (j = 0; j < ntypes; j++) {
|
for (i = 0; i < ntypes; i++) for (j = 0; j < ntypes; j++) {
|
||||||
fast_gamma_t* _noalias tab = &fast_gamma[i*ntypes*nr+j*nr];
|
fast_gamma_t* _noalias tab = &fast_gamma[i*ntypes*nr+j*nr];
|
||||||
if (type2rhor[i+1][j+1] >= 0) {
|
if (type2rhor[i+1][j+1] >= 0) {
|
||||||
|
|||||||
@ -221,9 +221,9 @@ void PairMultiLucy::allocate()
|
|||||||
memory->create(cutsq,nt,nt,"pair:cutsq");
|
memory->create(cutsq,nt,nt,"pair:cutsq");
|
||||||
memory->create(tabindex,nt,nt,"pair:tabindex");
|
memory->create(tabindex,nt,nt,"pair:tabindex");
|
||||||
|
|
||||||
memset(&setflag[0][0],0,nt*nt*sizeof(int));
|
memset(&setflag[0][0],0,sizeof(int)*nt*nt);
|
||||||
memset(&cutsq[0][0],0,nt*nt*sizeof(double));
|
memset(&cutsq[0][0],0,sizeof(double)*nt*nt);
|
||||||
memset(&tabindex[0][0],0,nt*nt*sizeof(int));
|
memset(&tabindex[0][0],0,sizeof(int)*nt*nt);
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
|
|||||||
@ -310,9 +310,9 @@ void PairMultiLucyRX::allocate()
|
|||||||
memory->create(cutsq,nt,nt,"pair:cutsq");
|
memory->create(cutsq,nt,nt,"pair:cutsq");
|
||||||
memory->create(tabindex,nt,nt,"pair:tabindex");
|
memory->create(tabindex,nt,nt,"pair:tabindex");
|
||||||
|
|
||||||
memset(&setflag[0][0],0,nt*nt*sizeof(int));
|
memset(&setflag[0][0],0,sizeof(int)*nt*nt);
|
||||||
memset(&cutsq[0][0],0,nt*nt*sizeof(double));
|
memset(&cutsq[0][0],0,sizeof(double)*nt*nt);
|
||||||
memset(&tabindex[0][0],0,nt*nt*sizeof(int));
|
memset(&tabindex[0][0],0,sizeof(int)*nt*nt);
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
|
|||||||
@ -919,7 +919,7 @@ void FixIMD::post_force(int /*vflag*/)
|
|||||||
|
|
||||||
if (imd_forces < length) { /* grow holding space for forces, if needed. */
|
if (imd_forces < length) { /* grow holding space for forces, if needed. */
|
||||||
memory->destroy(force_buf);
|
memory->destroy(force_buf);
|
||||||
force_buf = (void *) memory->smalloc(length*size_one,
|
force_buf = (void *) memory->smalloc((bigint)length*size_one,
|
||||||
"imd:force_buf");
|
"imd:force_buf");
|
||||||
}
|
}
|
||||||
imd_forces = length;
|
imd_forces = length;
|
||||||
@ -960,7 +960,7 @@ void FixIMD::post_force(int /*vflag*/)
|
|||||||
if (old_imd_forces < imd_forces) { /* grow holding space for forces, if needed. */
|
if (old_imd_forces < imd_forces) { /* grow holding space for forces, if needed. */
|
||||||
if (force_buf != nullptr)
|
if (force_buf != nullptr)
|
||||||
memory->sfree(force_buf);
|
memory->sfree(force_buf);
|
||||||
force_buf = memory->smalloc(imd_forces*size_one, "imd:force_buf");
|
force_buf = memory->smalloc((bigint)imd_forces*size_one, "imd:force_buf");
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
MPI_Bcast(force_buf, imd_forces*size_one, MPI_BYTE, 0, world);
|
MPI_Bcast(force_buf, imd_forces*size_one, MPI_BYTE, 0, world);
|
||||||
|
|||||||
@ -136,7 +136,7 @@ FixTTMMod::FixTTMMod(LAMMPS *lmp, int narg, char **arg) :
|
|||||||
gfactor1 = new double[atom->ntypes+1];
|
gfactor1 = new double[atom->ntypes+1];
|
||||||
gfactor2 = new double[atom->ntypes+1];
|
gfactor2 = new double[atom->ntypes+1];
|
||||||
// allocate 3d grid variables
|
// allocate 3d grid variables
|
||||||
total_nnodes = nxnodes*nynodes*nznodes;
|
total_nnodes = (bigint)nxnodes*nynodes*nznodes;
|
||||||
memory->create(nsum,nxnodes,nynodes,nznodes,"ttm/mod:nsum");
|
memory->create(nsum,nxnodes,nynodes,nznodes,"ttm/mod:nsum");
|
||||||
memory->create(nsum_all,nxnodes,nynodes,nznodes,"ttm/mod:nsum_all");
|
memory->create(nsum_all,nxnodes,nynodes,nznodes,"ttm/mod:nsum_all");
|
||||||
memory->create(sum_vsq,nxnodes,nynodes,nznodes,"ttm/mod:sum_vsq");
|
memory->create(sum_vsq,nxnodes,nynodes,nznodes,"ttm/mod:sum_vsq");
|
||||||
|
|||||||
@ -723,7 +723,7 @@ void FixPhonon::postprocess( )
|
|||||||
fwrite(&nucell, sizeof(int), 1, fp_bin);
|
fwrite(&nucell, sizeof(int), 1, fp_bin);
|
||||||
fwrite(&boltz, sizeof(double), 1, fp_bin);
|
fwrite(&boltz, sizeof(double), 1, fp_bin);
|
||||||
|
|
||||||
fwrite(Phi_all[0],sizeof(double),ntotal*fft_dim2*2,fp_bin);
|
fwrite(Phi_all[0],sizeof(double),(bigint)ntotal*fft_dim2*2,fp_bin);
|
||||||
|
|
||||||
fwrite(&TempAve, sizeof(double),1, fp_bin);
|
fwrite(&TempAve, sizeof(double),1, fp_bin);
|
||||||
fwrite(&basevec[0], sizeof(double),9, fp_bin);
|
fwrite(&basevec[0], sizeof(double),9, fp_bin);
|
||||||
|
|||||||
@ -302,10 +302,10 @@ int Allocate_Workspace( reax_system * /*system*/, control_params * control,
|
|||||||
|
|
||||||
// storage for reductions with multiple threads
|
// storage for reductions with multiple threads
|
||||||
#ifdef LMP_USER_OMP
|
#ifdef LMP_USER_OMP
|
||||||
workspace->CdDeltaReduction = (double *) scalloc(control->error_ptr, sizeof(double), total_cap*control->nthreads,
|
workspace->CdDeltaReduction = (double *) scalloc(control->error_ptr, sizeof(double), (rc_bigint)total_cap*control->nthreads,
|
||||||
"cddelta_reduce");
|
"cddelta_reduce");
|
||||||
|
|
||||||
workspace->forceReduction = (rvec *) scalloc(control->error_ptr, sizeof(rvec), total_cap*control->nthreads,
|
workspace->forceReduction = (rvec *) scalloc(control->error_ptr, sizeof(rvec), (rc_bigint)total_cap*control->nthreads,
|
||||||
"forceReduction");
|
"forceReduction");
|
||||||
|
|
||||||
workspace->valence_angle_atom_myoffset = (int *) scalloc(control->error_ptr, sizeof(int), total_cap, "valence_angle_atom_myoffset");
|
workspace->valence_angle_atom_myoffset = (int *) scalloc(control->error_ptr, sizeof(int), total_cap, "valence_angle_atom_myoffset");
|
||||||
|
|||||||
@ -113,7 +113,7 @@ void AngleHybrid::compute(int eflag, int vflag)
|
|||||||
|
|
||||||
const int nthreads = comm->nthreads;
|
const int nthreads = comm->nthreads;
|
||||||
if (comm->nthreads > 1) {
|
if (comm->nthreads > 1) {
|
||||||
const int nall = atom->nlocal + atom->nghost;
|
const bigint nall = atom->nlocal + atom->nghost;
|
||||||
if (eflag_atom)
|
if (eflag_atom)
|
||||||
memset(&eatom[0],0,nall*nthreads*sizeof(double));
|
memset(&eatom[0],0,nall*nthreads*sizeof(double));
|
||||||
if (vflag_atom)
|
if (vflag_atom)
|
||||||
|
|||||||
@ -113,7 +113,7 @@ void BondHybrid::compute(int eflag, int vflag)
|
|||||||
|
|
||||||
const int nthreads = comm->nthreads;
|
const int nthreads = comm->nthreads;
|
||||||
if (nthreads > 1) {
|
if (nthreads > 1) {
|
||||||
const int nall = atom->nlocal + atom->nghost;
|
const bigint nall = atom->nlocal + atom->nghost;
|
||||||
if (eflag_atom)
|
if (eflag_atom)
|
||||||
memset(&eatom[0],0,nall*nthreads*sizeof(double));
|
memset(&eatom[0],0,nall*nthreads*sizeof(double));
|
||||||
if (vflag_atom)
|
if (vflag_atom)
|
||||||
|
|||||||
@ -979,7 +979,7 @@ rendezvous_all2all(int n, char *inbuf, int insize, int inorder, int *procs,
|
|||||||
|
|
||||||
offsets[0] = 0;
|
offsets[0] = 0;
|
||||||
for (int i = 1; i < nprocs; i++)
|
for (int i = 1; i < nprocs; i++)
|
||||||
offsets[i] = offsets[i-1] + insize*procs_a2a[i-1];
|
offsets[i] = offsets[i-1] + (bigint)insize*procs_a2a[i-1];
|
||||||
|
|
||||||
bigint offset = 0;
|
bigint offset = 0;
|
||||||
for (int i = 0; i < n; i++) {
|
for (int i = 0; i < n; i++) {
|
||||||
@ -989,7 +989,8 @@ rendezvous_all2all(int n, char *inbuf, int insize, int inorder, int *procs,
|
|||||||
offset += insize;
|
offset += insize;
|
||||||
}
|
}
|
||||||
|
|
||||||
all2all1_bytes = nprocs*sizeof(int) + nprocs*sizeof(bigint) + n*insize;
|
all2all1_bytes = nprocs*sizeof(int) + nprocs*sizeof(bigint)
|
||||||
|
+ (bigint)n*insize;
|
||||||
|
|
||||||
} else {
|
} else {
|
||||||
procs_a2a = procs;
|
procs_a2a = procs;
|
||||||
@ -1085,7 +1086,7 @@ rendezvous_all2all(int n, char *inbuf, int insize, int inorder, int *procs,
|
|||||||
|
|
||||||
offsets[0] = 0;
|
offsets[0] = 0;
|
||||||
for (int i = 1; i < nprocs; i++)
|
for (int i = 1; i < nprocs; i++)
|
||||||
offsets[i] = offsets[i-1] + outsize*procs_a2a[i-1];
|
offsets[i] = offsets[i-1] + (bigint)outsize*procs_a2a[i-1];
|
||||||
|
|
||||||
bigint offset = 0;
|
bigint offset = 0;
|
||||||
for (int i = 0; i < nrvous_out; i++) {
|
for (int i = 0; i < nrvous_out; i++) {
|
||||||
@ -1096,7 +1097,7 @@ rendezvous_all2all(int n, char *inbuf, int insize, int inorder, int *procs,
|
|||||||
}
|
}
|
||||||
|
|
||||||
all2all2_bytes = nprocs*sizeof(int) + nprocs*sizeof(bigint) +
|
all2all2_bytes = nprocs*sizeof(int) + nprocs*sizeof(bigint) +
|
||||||
nrvous_out*outsize;
|
(bigint)nrvous_out*outsize;
|
||||||
|
|
||||||
} else {
|
} else {
|
||||||
procs_a2a = procs_rvous;
|
procs_a2a = procs_rvous;
|
||||||
|
|||||||
@ -258,7 +258,7 @@ void ComputeOrientOrderAtom::compute_peratom()
|
|||||||
|
|
||||||
double **x = atom->x;
|
double **x = atom->x;
|
||||||
int *mask = atom->mask;
|
int *mask = atom->mask;
|
||||||
memset(&qnarray[0][0],0,nmax*ncol*sizeof(double));
|
memset(&qnarray[0][0],0,sizeof(double)*nmax*ncol);
|
||||||
|
|
||||||
for (ii = 0; ii < inum; ii++) {
|
for (ii = 0; ii < inum; ii++) {
|
||||||
i = ilist[ii];
|
i = ilist[ii];
|
||||||
|
|||||||
10
src/memory.h
10
src/memory.h
@ -410,18 +410,18 @@ class Memory : protected Pointers {
|
|||||||
nbytes = ((bigint) sizeof(TYPE ***)) * n1;
|
nbytes = ((bigint) sizeof(TYPE ***)) * n1;
|
||||||
array = (TYPE ****) smalloc(nbytes,name);
|
array = (TYPE ****) smalloc(nbytes,name);
|
||||||
|
|
||||||
int i,j,k;
|
bigint i,j,k;
|
||||||
bigint m1,m2;
|
bigint m1,m2;
|
||||||
bigint n = 0;
|
bigint n = 0;
|
||||||
for (i = 0; i < n1; i++) {
|
for (i = 0; i < n1; i++) {
|
||||||
m2 = ((bigint) i) * n2;
|
m2 = i * n2;
|
||||||
array[i] = &plane[m2];
|
array[i] = &plane[m2];
|
||||||
for (j = 0; j < n2; j++) {
|
for (j = 0; j < n2; j++) {
|
||||||
m1 = ((bigint) i) * n2 + j;
|
m1 = i * n2 + j;
|
||||||
m2 = ((bigint) i) * n2*n3 + j*n3;
|
m2 = i * n2*n3 + j*n3;
|
||||||
plane[m1] = &cube[m2];
|
plane[m1] = &cube[m2];
|
||||||
for (k = 0; k < n3; k++) {
|
for (k = 0; k < n3; k++) {
|
||||||
m1 = ((bigint) i) * n2*n3 + j*n3 + k;
|
m1 = i * n2*n3 + j*n3 + k;
|
||||||
cube[m1] = &data[n];
|
cube[m1] = &data[n];
|
||||||
n += n4;
|
n += n4;
|
||||||
}
|
}
|
||||||
|
|||||||
Reference in New Issue
Block a user