ESPResSo 3.2.0-11-g9950804-git
Extensible Simulation Package for Soft Matter Research
debug.c
Go to the documentation of this file.
00001 /*
00002   Copyright (C) 2010,2011,2012,2013 The ESPResSo project
00003   Copyright (C) 2002,2003,2004,2005,2006,2007,2008,2009,2010 
00004     Max-Planck-Institute for Polymer Research, Theory Group
00005   
00006   This file is part of ESPResSo.
00007   
00008   ESPResSo is free software: you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation, either version 3 of the License, or
00011   (at your option) any later version.
00012   
00013   ESPResSo is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017   
00018   You should have received a copy of the GNU General Public License
00019   along with this program.  If not, see <http://www.gnu.org/licenses/>. 
00020 */
00021 /** \file debug.c
00022     Implements the malloc replacements as described in \ref debug.h "debug.h". 
00023 */
00024 
00025 #include <mpi.h>
00026 #include <signal.h>
00027 #include <stdio.h>
00028 #include <stdlib.h>
00029 #include <unistd.h>
00030 #include <string.h>
00031 #include "utils.h"
00032 #include "communication.h"
00033 #include "cells.h"
00034 #include "grid.h"
00035 #include "integrate.h"
00036 
00037 #if defined FORCE_CORE || defined MPI_CORE
00038 int regular_exit = 0;
00039 #else
00040 int regular_exit = 1;
00041 #endif
00042 static int core_done = 0;
00043 
00044 #ifdef ONEPART_DEBUG
00045 int check_id =  ONEPART_DEBUG_ID ;
00046 #endif
00047 
00048 #ifdef MEM_DEBUG
00049 
00050 #undef realloc
00051 #undef malloc
00052 #undef free
00053 
00054 void *__realloc(void *old, unsigned int size, char *where, int line)
00055 {
00056   void *ret;
00057   if (size == 0) {
00058     fprintf(stderr, "%d: free %p at %s:%d\n", this_node, old, where, line);
00059     free(old);
00060     return NULL;
00061   }
00062 
00063   ret = realloc(old, size);
00064   fprintf(stderr, "%d: realloc %p -> %p size %d at %s:%d\n", this_node, old, ret, size, where, line);
00065   return ret;
00066 }
00067 
00068 void *__malloc(unsigned int size, char *where, int line)
00069 {
00070   void *ret;
00071   if (size > 0)
00072     ret = malloc(size);
00073   else
00074     ret = NULL;
00075   fprintf(stderr, "%d: alloc %d -> %p at %s:%d\n", this_node, size, ret, where, line);
00076   return ret;
00077 }
00078 
00079 void __free(void *p, char *where, int line)
00080 {
00081   fprintf(stderr, "%d: free %p at %s:%d\n", this_node, p, where, line);
00082   free(p);
00083 }
00084 
00085 #endif
00086 
00087 void core()
00088 {
00089   if (!core_done && !regular_exit) {
00090     fprintf(stderr, "%d: forcing core dump on irregular exit (%d / %d) \n", this_node, core_done, regular_exit);
00091     /* this produced a compiler warning and got thrown out by GCC: */
00092     /* *(int *)0 = 0; */
00093     abort();
00094     /* doesn't work on AMD64 */
00095     /* kill(getpid(), SIGSEGV); */
00096     core_done = 1;
00097   }
00098 }
00099 
00100 void check_particle_consistency()
00101 {
00102   Particle *part;
00103   Cell *cell;
00104   int n, dir, c, p;
00105   int cell_part_cnt=0, ghost_part_cnt=0, local_part_cnt=0;
00106   int cell_err_cnt=0;
00107 
00108   /* checks: part_id, part_pos, local_particles id */
00109   for (c = 0; c < local_cells.n; c++) {
00110     cell = local_cells.cell[c];
00111     cell_part_cnt += cell->n;
00112     part = cell->part;
00113     for(int n=0; n<cell->n ; n++) {
00114       if(part[n].p.identity < 0 || part[n].p.identity > max_seen_particle) {
00115         fprintf(stderr,"%d: check_particle_consistency: ERROR: Cell %d Part %d has corrupted id=%d\n",
00116                 this_node,c,n,cell->part[n].p.identity);
00117         errexit();
00118       }
00119       for(dir=0;dir<3;dir++) {
00120         if(PERIODIC(dir) && (part[n].r.p[dir] < -ROUND_ERROR_PREC*box_l[dir] || part[n].r.p[dir] - box_l[dir] > ROUND_ERROR_PREC*box_l[dir])) {
00121           fprintf(stderr,"%d: check_particle_consistency: ERROR: illegal pos[%d]=%f of part %d id=%d in cell %d\n",
00122                   this_node,dir,part[n].r.p[dir],n,part[n].p.identity,c);
00123           errexit();
00124         }
00125       }
00126       if(local_particles[part[n].p.identity] != &part[n]) {
00127         fprintf(stderr,"%d: check_particle_consistency: ERROR: address mismatch for part id %d: local: %p cell: %p in cell %d\n",
00128                 this_node,part[n].p.identity,local_particles[part[n].p.identity],
00129                 &part[n],c);
00130         errexit();
00131         
00132       }
00133     }
00134   }
00135 
00136   for (c = 0; c < ghost_cells.n; c++) {
00137     cell = ghost_cells.cell[c];
00138     if(cell->n>0) {
00139       ghost_part_cnt += cell->n;
00140       fprintf(stderr,"%d: check_particle_consistency: WARNING: ghost_cell %d contains %d particles!\n",
00141               this_node,c,cell->n);
00142     }
00143   }
00144   CELL_TRACE(fprintf(stderr,"%d: check_particle_consistency: %d particles in cells, %d particles in ghost_cells.\n",
00145                      this_node,cell_part_cnt, ghost_part_cnt));
00146   /* checks: local particle id */
00147   for(n=0; n< max_seen_particle+1; n++) {
00148     if(local_particles[n] != NULL) {
00149       local_part_cnt ++;
00150       if(local_particles[n]->p.identity != n) {
00151         fprintf(stderr,"%d: check_particle_consistency: ERROR: local_particles part %d has corrupted id %d\n",
00152                 this_node,n,local_particles[n]->p.identity);
00153         errexit();
00154       }
00155     }
00156   }
00157   CELL_TRACE(fprintf(stderr,"%d: check_particle_consistency: %d particles in local_particles.\n",
00158                      this_node,local_part_cnt));
00159 
00160   /* EXIT on severe errors */
00161   if(cell_err_cnt>0) {
00162     fprintf(stderr,"%d: check_particle_consistency: %d ERRORS detected in cell structure!\n",this_node,cell_err_cnt);
00163     errexit();
00164   }
00165   if(local_part_cnt != cell_part_cnt) {
00166     fprintf(stderr,"%d: check_particle_consistency: ERROR: %d parts in cells but %d parts in local_particles\n",
00167             this_node,cell_part_cnt,local_part_cnt);
00168 
00169     for (c = 0; c < local_cells.n; c++) {
00170       for(p = 0; p < local_cells.cell[c]->n; p++)
00171         fprintf(stderr, "%d: got particle %d in cell %d\n", this_node, local_cells.cell[c]->part[p].p.identity, c);
00172     }
00173     
00174     for(p = 0; p < n_total_particles; p++)
00175       if (local_particles[p])
00176         fprintf(stderr, "%d: got particle %d in local_particles\n", this_node, p);
00177 
00178     if(ghost_part_cnt==0) errexit();
00179   }
00180   if(ghost_part_cnt>0) {
00181     fprintf(stderr,"%d: check_particle_consistency: ERROR: Found %d illegal ghost particles!\n",
00182             this_node,ghost_part_cnt);
00183     errexit();
00184   }
00185 }
00186 
00187 void check_particles()
00188 {
00189   Particle *part;
00190   int *is_here;
00191   Cell *cell;
00192   int n, dir, c, p;
00193   int cell_part_cnt=0, local_part_cnt=0;
00194   int cell_err_cnt=0;
00195   double skin2 = (skin != -1) ? skin/2 : 0;
00196 
00197   CELL_TRACE(fprintf(stderr, "%d: entering check_particles\n", this_node));
00198 
00199   /* check the consistency of particle_nodes */
00200   /* to this aim the array is broadcasted temporarily */
00201   if (this_node != 0)
00202     particle_node = malloc((max_seen_particle + 1)*sizeof(int));
00203   is_here = malloc((max_seen_particle + 1)*sizeof(int));
00204   memset(is_here, 0, (max_seen_particle + 1)*sizeof(int));
00205 
00206   MPI_Bcast(particle_node, max_seen_particle + 1, MPI_INT, 0, comm_cart);
00207 
00208   /* checks: part_id, part_pos, local_particles id */
00209   for (c = 0; c < local_cells.n; c++) {
00210     cell = local_cells.cell[c];
00211     cell_part_cnt += cell->n;
00212     part = cell->part;
00213     for(int n=0; n<cell->n ; n++) {
00214       if(part[n].p.identity < 0 || part[n].p.identity > max_seen_particle) {
00215         fprintf(stderr,"%d: check_particles: ERROR: Cell %d Part %d has corrupted id=%d\n",
00216                 this_node,c,n,cell->part[n].p.identity);
00217         errexit();
00218       }
00219 
00220       is_here[part[n].p.identity] = 1;
00221 
00222       for(dir=0;dir<3;dir++) {
00223         if(PERIODIC(dir) && (part[n].r.p[dir] < -skin2 || part[n].r.p[dir] > box_l[dir] + skin2)) {
00224           fprintf(stderr,"%d: check_particles: ERROR: illegal pos[%d]=%f of part %d id=%d in cell %d\n",
00225                   this_node,dir,part[n].r.p[dir],n,part[n].p.identity,c);
00226           errexit();
00227         }
00228       }
00229       if(local_particles[part[n].p.identity] != &part[n]) {
00230         fprintf(stderr,"%d: check_particles: ERROR: address mismatch for part id %d: local: %p cell: %p in cell %d\n",
00231                 this_node,part[n].p.identity,local_particles[part[n].p.identity],
00232                 &part[n],c);
00233         errexit();
00234       }
00235       if (particle_node[part[n].p.identity] != this_node) {
00236         fprintf(stderr,"%d: check_particles: ERROR: node for particle %d wrong\n",
00237                 this_node,part[n].p.identity);
00238         errexit();
00239       }
00240     }
00241   }
00242   CELL_TRACE(fprintf(stderr,"%d: check_particles: %d particles in local cells.\n",
00243                      this_node,cell_part_cnt));
00244 
00245   /* checks: local particle id */
00246   for(n = 0; n <= max_seen_particle; n++) {
00247     if(local_particles[n] != NULL) {
00248       local_part_cnt ++;
00249       if(local_particles[n]->p.identity != n) {
00250         fprintf(stderr,"%d: check_particles: ERROR: local_particles part %d has corrupted id %d\n",
00251                 this_node,n,local_particles[n]->p.identity);
00252         errexit();
00253       }
00254     }
00255   }
00256   CELL_TRACE(fprintf(stderr,"%d: check_particles: %d particles in local_particles.\n",
00257                      this_node,local_part_cnt));
00258 
00259   /* EXIT on severe errors */
00260   if(cell_err_cnt>0) {
00261     fprintf(stderr,"%d: check_particles: %d ERRORS detected in cell structure!\n",this_node,cell_err_cnt);
00262     errexit();
00263   }
00264 
00265   /* check whether the particles on my node are actually here */
00266   for (p = 0; p <= max_seen_particle; p++) {
00267     if (particle_node[p] == this_node) {
00268       if (!is_here[p]) {
00269         fprintf(stderr,"%d: check_particles: ERROR: particle %d on this node, but not in local cell\n", this_node, p);
00270       }
00271     }
00272   }
00273 
00274   free(is_here);
00275 
00276   if (this_node != 0) {
00277     free(particle_node);
00278     particle_node = NULL;
00279   }
00280   else {
00281     /* check whether the total count of particles is ok */
00282     c = 0;
00283     for (p = 0; p <= max_seen_particle; p++)
00284       if (particle_node[p] != -1) c++;
00285     if (c != n_total_particles) {
00286       fprintf(stderr,"%d: check_particles: #particles in particle_node inconsistent\n", this_node);
00287       errexit();
00288     }
00289     CELL_TRACE(fprintf(stderr,"%d: check_particles: %d particles in particle_node.\n",
00290                        this_node,c));
00291   }
00292   CELL_TRACE(fprintf(stderr, "%d: leaving check_particles\n", this_node));
00293 }
00294 
00295 void print_particle_positions()
00296 {
00297   int c,np,p;
00298   Cell *cell;
00299   Particle *part;
00300 
00301   for(c=0;c<n_cells;c++) {
00302     cell = &cells[c];
00303     part = cell->part;
00304     np = cell->n;
00305     if(np>0) {
00306       fprintf(stderr,"%d: cell %d contains %d particles:\n",this_node,c,np);
00307       for(p=0;p<np;p++) {
00308         fprintf(stderr,"%d: c%d p%d (id%d) pos %f %f %f\n",this_node,c,p,part[p].p.identity,part[p].r.p[0],part[p].r.p[1],part[p].r.p[2]);
00309       }
00310     }
00311   }    
00312 }
00313 
00314 void print_particle_forces()
00315 {
00316   int c,np,p;
00317   Cell *cell;
00318   Particle *part;
00319 
00320   for(c=0;c<n_cells;c++) {
00321     cell = &cells[c];
00322     part = cell->part;
00323     np = cell->n;
00324     if(np>0) {
00325       fprintf(stderr,"%d: cell %d contains %d particles:\n",this_node,c,np);
00326       for(p=0;p<np;p++) {
00327         fprintf(stderr,"%d: c%d p%d (id%d) force %f %f %f\n",this_node,c,p,part[p].p.identity,part[p].f.f[0],part[p].f.f[1],part[p].f.f[2]);
00328       }
00329     }
00330   }    
00331 }