Timothy  0.9
Tissue Modelling Framework
 All Data Structures Files Functions Variables Typedefs Macros
stats.c
Go to the documentation of this file.
1 /* **************************************************************************
2  * This file is part of Timothy
3  *
4  * Copyright (c) 2014/15 Maciej Cytowski
5  * Copyright (c) 2014/15 ICM, University of Warsaw, Poland
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
20  *
21  * *************************************************************************/
22 
23 #include <stdio.h>
24 #include <stdlib.h>
25 #include <mpi.h>
26 #include <math.h>
27 #include <inttypes.h>
28 
29 #include "global.h"
30 
39 {
40  int c;
41  double sum = 0.0, mean;
42  double globalMaxDens;
43  double globalMinDens;
44  double deviation;
45 
46  statistics.maxdens = 0;
47  statistics.mindens = 1024;
48 
49  for (c = 0; c < lnc; c++) {
50  sum += cells[c].density;
52  (cells[c].density >
55  (cells[c].density <
57  }
58 
59  MPI_Allreduce(&sum, &mean, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
60  MPI_Allreduce(&statistics.maxdens, &globalMaxDens, 1, MPI_DOUBLE,
61  MPI_MAX, MPI_COMM_WORLD);
62  MPI_Allreduce(&statistics.mindens, &globalMinDens, 1, MPI_DOUBLE,
63  MPI_MIN, MPI_COMM_WORLD);
64  statistics.maxdens = globalMaxDens;
65  statistics.mindens = globalMinDens;
66 
67  mean /= nc;
68  statistics.densavg = mean;
69  sum = 0.0;
70 
71  for (c = 0; c < lnc; c++)
72  sum += (cells[c].density - mean) * (cells[c].density - mean);
73 
74  MPI_Allreduce(&sum, &(statistics.densdev), 1, MPI_DOUBLE, MPI_SUM,
75  MPI_COMM_WORLD);
76 
79 
80  if (MPIrank == 0)
81  printf("%12s%10.4lf%10.4lf%10.4lf%10.4lf\n", "Density ",
82  globalMinDens, globalMaxDens, mean, statistics.densdev);
83 }
84 
89 {
90  double globalMinDist;
91  double globalMaxDist;
92 
93  MPI_Allreduce(&statistics.mindist, &globalMinDist, 1, MPI_DOUBLE,
94  MPI_MIN, MPI_COMM_WORLD);
95 
96  statistics.mindist = globalMinDist;
97 
98  if (MPIrank == 0)
99  printf("%12s%10.4lf%10s%10s%10s\n", "Distance ", globalMinDist, "-",
100  "-", "-");
101 }
102 
107 {
108  MPI_Reduce(&statistics.minvel, &globalMinVel, 1, MPI_DOUBLE, MPI_MIN, 0,
109  MPI_COMM_WORLD);
110  MPI_Reduce(&statistics.maxvel, &globalMaxVel, 1, MPI_DOUBLE, MPI_MAX, 0,
111  MPI_COMM_WORLD);
112  if (MPIrank == 0)
113  printf("%12s%10.4lf%10.4lf%10s%10s\n", "Velocity ", globalMinVel,
114  globalMaxVel, "-", "-");
115 }
116 
121 {
122  double globalMaxSize, globalMinSize;
123 
124  MPI_Reduce(&statistics.maxsize, &globalMaxSize, 1, MPI_DOUBLE, MPI_MAX,
125  0, MPI_COMM_WORLD);
126  MPI_Reduce(&statistics.minsize, &globalMinSize, 1, MPI_DOUBLE, MPI_MIN,
127  0, MPI_COMM_WORLD);
128 
129  if (MPIrank == 0)
130  printf("%12s%10.4lf%10.4lf%10s%10s\n", "Size ", globalMinSize,
131  globalMaxSize, "-", "-");
132 }
133 
138 {
139  if (MPIrank == 0) {
140  printf("%12s%12s%12s%12s%12s%12s\n", "Cell phase ", "G0", "G1", "S",
141  "G2", "M");
142  printf("%12s%12" PRId64 "%12" PRId64 "%12" PRId64 "%12" PRId64 "%12"
143  PRId64 "\n", "N. of cells", g0nc, g1nc, snc, g2nc, mnc);
144  printf("\n%16s%12" PRId64 "\n", "Healthy cells ", nc - cnc - nnc);
145  printf("%16s%12" PRId64 "\n", "Cancer cells ", cnc);
146  printf("%16s%12" PRId64 "\n", "Necrotic cells ", nnc);
147  }
148 }
149 
154 {
155  if (MPIrank == 0)
156  printf("%12s%10s%10s%10s%10s\n", "", "Min", "Max", "Avg", "Dev");
159  /*statisticsVelocity(); */
160  statisticsSize();
161  if (MPIrank == 0)
162  printf("\n");
164 }
double density
Definition: global.h:75
int MPIrank
Definition: global.h:134
#define lnc
Definition: global.h:102
double densavg
Definition: global.h:265
double densdev
Definition: global.h:264
double maxdens
Definition: global.h:262
double mindist
Definition: global.h:259
double maxsize
Definition: global.h:258
contains the most important global variables, arrays and defines
void statisticsSize()
Definition: stats.c:120
double minsize
Definition: global.h:257
void statisticsVelocity()
Definition: stats.c:106
void statisticsDistance()
Definition: stats.c:88
struct statisticsData statistics
Definition: global.h:268
#define g1nc
Definition: global.h:95
void statisticsPhases()
Definition: stats.c:137
#define snc
Definition: global.h:96
struct cellData * cells
Definition: global.h:82
void statisticsDensity()
Definition: stats.c:38
double maxvel
Definition: global.h:260
void statisticsPrint()
Definition: stats.c:153
#define g0nc
Definition: global.h:94
#define nnc
Definition: global.h:100
#define cnc
Definition: global.h:99
#define nc
Definition: global.h:93
double globalMinVel
Definition: global.h:270
#define g2nc
Definition: global.h:97
#define mnc
Definition: global.h:98
double globalMaxVel
Definition: global.h:271
double mindens
Definition: global.h:263
double minvel
Definition: global.h:261