#include <fstream>
#include <sstream>
#include <cstring>
#include <functional>
#include <getopt.h>
#include <cerrno>
#include "mesh.h"
#include "parameter.h"
#include "random.h"
#include "pi.h"
#include "cellitem.h"
#ifdef QTGRAPHICS
#include "canvas.h"
#include "cell.h"
#include <qwidget.h>
#include <q3process.h>
#include <qapplication.h>
#include <QDesktopWidget>
#include <QGraphicsScene>
#include <QMessageBox>
//Added by qt3to4:
#include <QMouseEvent>
#endif
#include <unistd.h>
#include <q3textstream.h>
#ifdef USE_LICENSE
#include "license.h"
#endif
#ifdef HAVE_QWT
#include "data_plot.h"
#endif
#include <QPalette>
#include <QBrush>
#include <QToolTip>
//#include "reactions.h"
//#include "reactions_auxacc.h"
//#include "reactions_pinaux3.h"
//#include "reactions_aux1.h"
//#define REACTIONS_HEADER "reactions_growth.h"
//#include "reactions_tips_nor.h"
#define _xstr_(s) _str_(s)
#define _str_(s) #s
#include _xstr_(REACTIONS_HEADER)
extern Parameter par;
MainBase *main_window = 0;
double auxin_account = 0.;
#ifdef XFIGGRAPHICS
#define TIMESTEP double Graphics::TimeStep(void)
#endif
//ofstream debug_stream("debug.log");
class PrintNode {
public:
void operator() (const Node &n) const
{
cerr << n.Index() << ": " << n << endl;
}
};
class EdgeSource {
public:
void operator() (Cell &c) {
if (c.AtBoundaryP()) {
cerr << "Cell " << c.Index() << " is a source cell.\n";
c.SetSource(0,par.source);
} else {
cerr << "Cell " << c.Index() << " is _not_ a source cell.\n";
}
}
};
/*class SetArea {
public:
void operator() (Cell &c) const {
c.SetTargetArea(4*c.CalcArea());
c.SetTargetLength(4*c.Length());
}
};*/
#ifdef XFIGGRAPHICS
class XFigCell {
public:
void operator() (Cell &c,std::ostream &os) const {
c.XFigPrint(os);
}
};
#endif
class CellInfo {
public:
void operator() (Cell &c,std::ostream &os) const {
os << "Cell " << c.index << " says: " << endl;
os << "c.nodes.size() = " << c.nodes.size() << endl;
for (list<Node *>::iterator i=c.nodes.begin();
i!=c.nodes.end();
i++) {
cerr << (*i)->Index() << " ";
}
cerr << endl;
}
};
double PINSum(Cell &c) {
return c.Chemical(1) + c.SumTransporters(1);// + c.ReduceCellAndWalls<double>( complex_PijAj );
}
#ifdef QTGRAPHICS
class DrawCell {
public:
/*void operator() (Cell &c,Graphics &g) const {
c.Draw(g);
}*/
void operator() (Cell &c,QGraphicsScene &canvas, MainBase &m) const {
if (m.ShowBorderCellsP() || c.Boundary()==Cell::None) {
if (!m.ShowBoundaryOnlyP() && !m.HideCellsP())
if (m.ShowToolTipsP()) {
QString info_string=QString("Cell %1, chemicals: ( %2, %3, %4, %5, %6)\n %7 of PIN1 at walls.\n Area is %8\n PIN sum is %9\n Circumference is %10\n Boundary type is %11").arg(c.Index()).arg(c.Chemical(0)).arg(c.Chemical(1)).arg(c.Chemical(2)).arg(c.Chemical(3)).arg(c.Chemical(4)).arg(c.SumTransporters(1)).arg(c.Area()).arg(PINSum(c)).arg(c.Circumference()).arg(c.BoundaryStr());
info_string += "\n" + c.printednodelist();
c.Draw(&canvas, info_string);
} else {
c.Draw(&canvas);
}
if (m.ShowCentersP())
c.DrawCenter(&canvas);
// if (m.ShowMeshP())
// c.DrawNodes(&canvas);
if (m.ShowFluxesP())
c.DrawFluxes(&canvas, par.arrowsize);
//c.DrawTriangles(canvas);
}
//cerr << "Area of cell " << c.Index() << " is " << c.Area() << endl;
}
};
#endif
const char *xfig_header = "#FIG 3.2\nLandscape\nCenter\nInches\nLetter\n100.00\nSingle\n-2\n1200 2\n";
Mesh mesh;
bool batch=false;
#ifdef QTGRAPHICS
void MainBase::Plot(int resize_stride) {
clear();
static int count=0;
if (resize_stride) {
if ( !((++count)%resize_stride) ) {
FitLeafToCanvas();
}
}
mesh.LoopCells(DrawCell(),canvas,*this);
if (ShowNodeNumbersP())
mesh.LoopNodes( bind2nd (mem_fun_ref ( &Node::DrawIndex), &canvas ) ) ;
if (ShowCellNumbersP())
mesh.LoopCells( bind2nd (mem_fun_ref ( &Cell::DrawIndex), &canvas ) ) ;
if (ShowCellAxesP())
mesh.LoopCells( bind2nd (mem_fun_ref ( &Cell::DrawAxis), &canvas ) );
if (ShowCellStrainP())
mesh.LoopCells( bind2nd (mem_fun_ref ( &Cell::DrawStrain), &canvas ) );
if (ShowWallsP())
//mesh.LoopCells( bind2nd (mem_fun_ref( &Cell::DrawWalls), &canvas ) );
mesh.LoopWalls( bind2nd( mem_fun_ref( &Wall::Draw ), &canvas ) );
if (ShowApoplastsP())
mesh.LoopWalls( bind2nd( mem_fun_ref( &Wall::DrawApoplast ), &canvas ) );
if (ShowMeshP())
mesh.DrawNodes(&canvas);
if (ShowBoundaryOnlyP())
mesh.DrawBoundary(&canvas);
// check to see if chemicals emitted by leaf edge diffuse evenly. Yes, they do.
/*
Vector mesh_centroid = mesh.Centroid();
Vector first_moment = mesh.FirstConcMoment(0);
QCanvasEllipse *disk1 = 0;
if (disk1==0) {
disk1 = new QCanvasEllipse ( 50, 50, &canvas );
disk1->setBrush( QColor("red") );
disk1->setZ(5);
disk1->show();
}
QCanvasEllipse *disk2 = 0;
if (disk2==0) {
disk2=new QCanvasEllipse ( 50, 50, &canvas );
disk2->setBrush( QColor("blue") );
disk2->setZ(7);
disk2->show();
}
Vector offset = mesh.Offset();
double factor = mesh.Factor();
disk1 -> move((offset.x+mesh_centroid.x)*factor,(offset.y+mesh_centroid.y)*factor);
disk2 -> move((offset.x+first_moment.x)*factor,(offset.y+first_moment.y)*factor);
*/
if ( ( batch || MovieFramesP() )) {
static int frame = 0;
// frame numbers are sequential for the most frequently written file type.
// for the less frequently written file type they match the other type
if (!(count%par.storage_stride) ) {
stringstream fname;
fname << par.datadir << "/leaf.";
fname.fill('0');
fname.width(6);
//fname << mesh.getTime() << ".xml";
fname << frame << ".jpg";
if (par.storage_stride <= par.xml_storage_stride) {
frame++;
}
// Write high-res JPG snapshot every plot step
Save(fname.str().c_str(), "JPEG",1024,768);
}
if (!(count%par.xml_storage_stride)) {
stringstream fname;
fname << par.datadir << "/leaf.";
fname.fill('0');
fname.width(6);
fname << frame << ".xml";
if (par.xml_storage_stride < par.storage_stride) {
frame++;
}
// Write XML file every ten plot steps
mesh.XMLSave(fname.str().c_str(), XMLSettingsTree());
}
}
}
#endif
//double areaunit=0;
void Cell::Flux(double *flux, double *D) {
// Algorithm according to Rudge & Haseloff 2005
// (will we need to take cell area into account?)
// For the time being, we don't: assume cell area is
// mainly determined by vacuole.
// loop over cell edges
for (int c=0;c<Cell::nchem;c++) flux[c]=0.;
for (list<Wall *>::iterator i=walls.begin();
i!=walls.end();
i++) {
// leaf cannot take up chemicals from environment ("no flux boundary")
if ((*i)->c2->BoundaryPolP()) continue;
// calculate edge length
// (will later be updated during node displacement for efficiency)
//double edge_length = (m->nodes[(*i)->n1]-m->nodes[(*i)->n2]).Norm();
// flux depends on edge length and concentration difference
for (int c=0;c<Cell::nchem;c++) {
double phi = (*i)->length * ( D[c] ) * ( (*i)->c2->chem[c] - chem[c] );
if ((*i)->c1!=this) {
cerr << "Warning, bad cells boundary: " << (*i)->c1->index << ", " << index << endl;
}
flux[c] += phi;
}
}
}
INIT {
//Cell &circle=mesh.CircularCell(0,0,10,10);
if (leaffile) {
xmlNode *settings;
mesh.XMLRead(leaffile, &settings);
main_window->XMLReadSettings(settings);
xmlFree(settings);
main_window->UserMessage(QString("Ready. Time is %1").arg(mesh.getTimeHours().c_str()));
} else {
//Cell &circle=mesh.LeafPrimordium(10,50);
Cell &circle=mesh.CircularCell(0,0,10,10);
//circle.SetChemical(1,0);
// petiole is auxin sink
//mesh.getCell(1).SetSource(0,0);
// petiole becomes provascular cell
//mesh.getCell(1).SetChemical(1,0.6);
circle.SetTargetArea(circle.CalcArea());
//circle.SetChemical(2,1e-60);
mesh.SetBaseArea();
circle.SetChemical(1, par.Pi_tot );
circle.SetChemical(0, 0.);
//mesh.LoopCells(EdgeSource());
}
}
TIMESTEP {
static int i=0;
static int t=0;
static int ncells;
if (!batch) {
UserMessage(QString("Time: %1").arg(mesh.getTimeHours().c_str()),0);
}
//if (!(i%1000))
//Save(QString("leaf%1.eps").arg(i).ascii(), "PDF");
/* static ofstream *kinematic=0;
if (kinematic == 0) {
stringstream kinematic_fname;
kinematic_fname << par.datadir << "/kinematic.dat";
cerr << "Writing kinematic data to " << kinematic_fname.str() << endl;
kinematic = new ofstream(kinematic_fname.str().c_str());
}
*kinematic << t << " " << i << " " << mesh.Area() << " " << mesh.NCells() << endl;
*/
ncells=mesh.NCells();
double dh;
//const double en_threshold = 1;
//mesh.RandomlyLoopCells( mem_fun_ref(&Cell::Displace ));
//static bool dumpflag=false;
if(DynamicCellsP()) {
dh = mesh.DisplaceNodes();
//static ofstream enfile("energy.dat");
//enfile << i << " " << dh << "\n";
// Only allow for node insertion, cell division and cell growth
// if the system has equillibrized
// i.e. cell wall tension equillibrization is much faster
// than biological processes, including division, cell wall yielding
// and cell expansion
mesh.InsertNodes(); // (this amounts to cell wall yielding)
if ( (-dh) < par.energy_threshold) {
mesh.IncreaseCellCapacityIfNecessary();
mesh.LoopCurrentCells(CellHouseKeeping()); // this includes cell division
// Reaction diffusion
TransportFunction *transport_f = new AuxinTransport();
CellReaction *cellreaction_f = new AuxinAndDifferentiation();
WallReaction *wall_f = new PIN1Localization();
mesh.ReactDiffuse_New(transport_f, cellreaction_f, wall_f, par.rd_dt);
t++;
Plot();
}
} else {
TransportFunction *transport_f = new AuxinTransport();
CellReaction *cellreaction_f = new AuxinAndDifferentiation();
WallReaction *wall_f = new PIN1Localization();
mesh.ReactDiffuse_New(transport_f, cellreaction_f, wall_f, par.rd_dt);
Plot();
}
i++;
return mesh.getTime();
}
void Main::Divide(void) {
static Vector axis(1,0,0);
mesh.IncreaseCellCapacityIfNecessary();
mesh.LoopCurrentCells( bind2nd( mem_fun_ref(&Cell::DivideOverAxis),axis) );
axis=axis.Perp2D();
//mesh.LoopCells( mem_fun_ref(&Cell::Divide) );
Plot();
}
/* Called if a cell is clicked */
void Cell::OnClick(QMouseEvent *e) {
cerr << "Calling OnClick()\n";
/* cerr << *m->boundary_polygon;
//m->RepairBoundaryPolygon();
//SetChemical(0,par.cellsource);
cerr << "Mesh's centroid: " << m->Centroid() << endl;
cerr << "First moment of chem[0]: " << m->FirstConcMoment(0) << endl;
*/
/*
cerr << "Cell: " << Index() << ", has walls: ";
for (list<Wall *>::const_iterator w=walls.begin();
w!=walls.end();
w++) {
cerr << **w << " ";
}
cerr << endl;
*/
/* cerr << "Cell: " << Index() << ", has neighbors: ";
for (list<Cell *>::const_iterator c=neighbors.begin();
c!=neighbors.end();
c++) {
cerr << (*c)->Index() << " ";
}
cerr << endl;*/
/*
mesh.PrintWallList();
*/
double circ=0.;
for (list<Wall *>::const_iterator w=walls.begin();
w!=walls.end();
w++) {
cerr << (*w)->N1()->Index() << "->" << (*w)->N2()->Index() << " = " << (*w)->Length() << endl;
circ += (*w)->Length();
}
cerr << "Circ is " << circ << endl;
if (e->button() == Qt::MidButton) {
double sum_pin_bef = getMesh().CalcProtCellsWalls(1);
QString message(QString("Dividing Cell %1").arg(index));
((Main *)main_window)->UserMessage(message);
cerr << message.toStdString();
Divide();
((Main *)main_window)->Plot();
double sum_pin_aft = getMesh().CalcProtCellsWalls(1);
cerr << "Sum PIN1, before: " << sum_pin_bef << ", after: " << sum_pin_aft << endl;
return;
}
if (e->button() == Qt::LeftButton) {
//double sum_walls=0.;
cerr << "Wall lengths: ";
for (list<Wall *>::const_iterator w=walls.begin();
w!=walls.end();
w++) {
//cerr << (*w)->getTransporter(this, 1) << " ";
cerr << (*w)->Length() << " ";
//sum_walls+=(*w)->getTransporter(this, 1); //* (*w)->Length();
}
cerr << ", Chemical(1) = " << Chemical(1) << ", sum = " << SumTransporters(1) + Chemical(1) << endl;
QString message;
message=QString("Cell %1 has chemicals ( %2, %3, %4, %5), and it has %6 of PIN1 at its walls. Area is %7").arg(Index()).arg(chem[0]).arg(chem[1]).arg(chem[2]).arg(chem[3]).arg(SumTransporters(1)).arg(Area());
((Main *)main_window)->UserMessage(message);
/* cerr << "Cell " << index << " has chemicals ( " << chem[0] << ", " << chem[1] << " ) " << endl;
cerr << "Cell " << index << "'s amount of PIN1 at the walls: " << SumTransporters(1) << endl;*/
#ifdef QTGRAPHICS
/* double sum_PIN1 = Chemical(4);
for (list<Wall *>::const_iterator i = walls.begin();
i!=walls.end();
i++) {
sum_PIN1 += (*i)->getTransporter(this, 0);
}
*/
//main_window->UserMessage(QString("Concentration chemical 0 of cell %1 = %2").arg(Index()).arg(chem[0]));
//main_window->UserMessage(QString("Target area of cell %1 = %2").arg(Index()).arg(TargetArea()));
//main_window->UserMessage(QString("Sum PIN1 of cell %1 = %2 (intracellular %3)").arg(Index()).arg(sum_PIN1).arg(Chemical(4)));
/* QString message;
message=QString("Cell %1's nodes are:").arg(Index());
for (list<int>::iterator it=neighbors.begin();
it!=neighbors.end();
it++) {
message += QString(" %2").arg(*it);
}
main_window->UserMessage(message);*/
//SetWallLengths();
//main_window->UserMessage(QString("Cell %1 apoptoses.").arg(Index()));
//Apoptose();
//m->CleanUpCellNodeLists();
#ifdef HAVE_QWT
QStringList curvenames;
curvenames += QString("Auxin");
curvenames += QString("PIN1");
curvenames += QString("auxill. prot");
curvenames += QString("Wall PIN1");
PlotDialog *plot = new PlotDialog((Main *)main_window, QString("Monitor for Cell %1").arg(Index()), curvenames);
QObject::connect(this, SIGNAL(ChemMonValue(double, double *)),
plot, SLOT(AddValue(double,double *)));
#endif
} else {
//Divide();
}
#endif
}
void Wall::OnWallInsert(void) {
// NOTE: THIS FUNCTION IS CALLED AFTER Cell::OnDivide();
// After division, we put some PIN1 on the walls, to prevent quick
// redistribution of PIN1s and consequent auxin accumulation after
// division.
// make sure dPij/dt = 0 at both sides of the new wall.
// True for:
// Pij = k1/k2 * A_j * L_ij * P_i
// Pji = k1/k2 * A_i * L_ij * P_j
//transporters1[1] = (par.k1/par.k2) * c2->Chemical(0) * length * c1->Chemical(1);
//transporters2[1] = (par.k1/par.k2) * c1->Chemical(0) * length * c2->Chemical(1);
//transporters1[1]=transporters2[1]=0.;
// cerr << "Length of new wall is " << length << endl;
//cerr << "PIN1 is [ " << transporters1[1] << " | " << transporters2[1] << "] ";
}
int main(int argc,char **argv) {
try {
//if (argc<2 || strstr(argv[1],".par")==0) {
// throw "Usage: Leaf [parfile].par";
//}
//par.Read(argv[1]);
// parse command-line options
int c;
//int digit_optind = 0;
char *parfile=0;
char *leaffile=0;
/* ofstream args("args.txt");
for (int i=0;i<argc;i++) {
args << argv[i] << endl;
}*/
while (1) {
//int this_option_optind = optind ? optind : 1;
int option_index = 0;
static struct option long_options[] = {
{"batch", 0, 0, 0},
{"par", 1, 0, 0},
{"leaffile", 2, 0, 0}
};
// short option 'p' creates trouble for non-commandline usage on MacOSX. Option -p changed to -P (capital)
static char *short_options = "bP:l";
c = getopt_long (argc, argv, "bP:l:",
long_options, &option_index);
if (c == -1)
break;
if (c==0) {
printf ("option %s", long_options[option_index].name);
if (optarg)
printf (" with arg %s", optarg);
printf ("\n");
c = short_options[option_index];
}
switch (c) {
case 'b':
cerr << "Running in batch mode\n";
batch=true;
break;
case 'p':
parfile=strdup(optarg);
if (!parfile) {
throw("Out of memory");
}
printf ("parameter file is '%s'\n", parfile);
break;
case 'l':
leaffile=strdup(optarg);
if (!leaffile) {
throw("Out of memory");
}
printf("Reading leaf state file '%s'\n", leaffile);
break;
case '?':
break;
default:
printf ("?? getopt returned character code 0%o ??\n", c);
}
}
if (optind < argc) {
printf ("non-option ARGV-elements: ");
while (optind < argc)
printf ("%s ", argv[optind++]);
printf ("\n");
}
#ifdef X11GRAPHICS
X11Graphics g(par.sizex,par.sizey);
g.ChangeTitle("Virtual Leaf");
#endif
#ifdef QTGRAPHICS
//QApplication app(argc, argv);
//argc=1;
QCoreApplication *app;
QGraphicsScene canvas(0,0,8000,6000);
QPalette tooltippalette = QToolTip::palette();
QColor transparentcolor = QColor(tooltippalette.brush(QPalette::Window).color());
//transparentcolor.setAlphaF(0.9);
tooltippalette.setBrush (QPalette::Window, QBrush (transparentcolor) );
QToolTip::setPalette( tooltippalette );
if (batch) {
// Note: QCoreApplication allows for command line applications, independent of X-Server on Unix.
// Allows for running on cluster
app = new QCoreApplication(argc,argv);
main_window=new MainBase(canvas, mesh);
} else {
app = new QApplication(argc,argv);
main_window=new Main(canvas, mesh);
((Main *)main_window)->resize( ((Main *)main_window)->sizeHint());
}
if (!batch) {
if ( QApplication::desktop()->width() > ((Main *)main_window)->width() + 10
&& QApplication::desktop()->height() > ((Main *)main_window)->height() +30 )
((Main *)main_window)->show();
else
((Main *)main_window)->showMaximized();
}
#ifdef USE_LICENSE
CheckLicenseFile();
#endif
canvas.setSceneRect(QRectF());
if (!batch) {
QObject::connect( qApp, SIGNAL(lastWindowClosed()), qApp, SLOT(quit()) );
}
// return app.exec();
#endif
if (parfile)
par.Read(parfile);
//else
// par.Read("default.par");
Seed(par.rseed);
#ifdef XFIGGRAPHICS
Graphics g;
#endif
main_window->Init(leaffile);
// if we have supplied a leaffile AND a parfile, reread the parameter file,
// because the leaffile also contains parameter information
if (leaffile && parfile) {
par.Read(parfile);
}
Cell::SetMagnification(10);
Cell::setOffset(0,0);
main_window->FitLeafToCanvas();
#ifdef QTGRAPHICS
// check if leaf crosses boundary
/* Vector bbll,bbur;
mesh.BoundingBox(bbll,bbur);
double scale_x = canvas.width()/(bbur.x-bbll.x);
double scale_y = canvas.height()/(bbur.y-bbll.y);
Cell::Scale(scale_x<scale_y ? scale_x:scale_y); // smallest of scale_x and scale_y
double offset_x = (canvas.width()/Cell::Magnification()-(bbur.x-bbll.x))/2.;
double offset_y = (canvas.height()/Cell::Magnification()-(bbur.y-bbll.y))/2.;
Cell::setOffset(offset_x, offset_y);
*/
/* Vector bbll,bbur;
mesh.BoundingBox(bbll,bbur);
cerr << "bbll = " << bbll << endl;
cerr << "bbur = " << bbur << endl;
double cw = (bbur.x - bbll.x);
double ch = (bbur.y - bbll.y);
double factor = canvas.width()/(2*cw);
cerr << "factor = " << factor << ", width = " << canvas.width() << endl;
cerr << "Size is " << cw << " x " << ch << endl;
canvas.resize((int)(2*cw*factor),(int)(2*ch*factor));
Cell::SetMagnification(factor);
//main_window->scale(factor);
//Cell::Translate(-bbll.x+cw/2,-bbll.y+ch/2);
cerr << -bbll.x+cw/2 << ", " << -bbll.y+ch/2 << ", " << canvas.width() << ", " << canvas.height() << endl;
Cell::Translate((-bbll.x+cw/2),(-bbll.y+ch/2));
*/
main_window->Plot();
#endif
#ifdef QTGRAPHICS
if (batch) {
//main_window->toggleMovieFrames();
//main_window->toggleShowNodes();
//for (int i=0;i<par.nit;i++) {
double t=0.;
do {
t = main_window->TimeStep();
/* cerr << endl << "***********************************************" << endl;
cerr << "Time is " << t << endl;
cerr << "par.maxt = " << par.maxt << endl;
cerr << endl << "***********************************************" << endl;*/
} while (t < par.maxt);
} else
return app->exec();
#else
//for (int i=0;i<par.nit;i++) {
do {
t= g.TimeStep();
} while (t < par.maxt);
#endif
} catch (const char *message) {
if (batch) {
cerr << "Exception caught:" << endl;
cerr << message << endl;
abort();
} else {
QString qmess=QString("Exception caught: %1").arg(message);
QMessageBox::critical(0, "Critical Error", qmess);
abort();
}
} catch (ios_base::failure) {
stringstream error_message;
error_message << "I/O failure: " << strerror(errno);
if (batch) {
cerr << error_message.str() <<endl;
abort();
} else {
QString qmess(error_message.str().c_str());
QMessageBox::critical(0, "I/O Error", qmess );
abort();
}
}
}
extern double auxin_account;
const int Cell::nchem = 4;
class AuxinTransport : public TransportFunction {
public:
virtual void operator()(Wall *w, double *dchem_c1, double *dchem_c2) {
// leaf edge is const source of auxin
// (Neumann boundary condition: we specify the influx)
if (w->C2()->BoundaryPolP()) {
if (w->AuxinSource()) {
double aux_flux = par.leaf_tip_source * w->Length();
dchem_c1[0]+= aux_flux;
// dchem_c2 is undefined..!
return;
} else {
if (w->AuxinSink()) {
// efflux into Shoot Apical meristem
// we assume all PINs are directed towards shoot apical meristem
dchem_c1[0] -= par.sam_efflux * w->C1()->Chemical(0) / (par.ka + w->C1()->Chemical(0));
return;
} else
return;
}
}
if (w->C1()->BoundaryPolP()) {
if (w->AuxinSource()) {
double aux_flux = par.leaf_tip_source * w->Length();
dchem_c2[0] += aux_flux;
// dchem_c1 is undefined...!
return;
} else {
if (w->AuxinSink()) {
// efflux into Shoot Apical meristem
// we assume all PINs are directed towards shoot apical meristem
dchem_c2[0] -= par.sam_efflux * w->C2()->Chemical(0) / (par.ka + w->C2()->Chemical(0));
return;
} else
return;
}
}
// Passive fluxes (Fick's law)
// only auxin flux now
// flux depends on edge length and concentration difference
int c=0;
double phi = w->Length() * ( par.D[c] ) * ( w->C2()->Chemical(c) - w->C1()->Chemical(c) );
dchem_c1[c] += phi;
dchem_c2[c] -= phi;
// Active fluxes (PIN1 mediated transport)
// (Transporters measured in moles, here)
// efflux from cell 1 to cell 2
double trans12 = ( par.transport * w->Transporters1(1) * w->C1()->Chemical(0) / (par.ka + w->C1()->Chemical(0)) );
// efflux from cell 2 to cell 1
double trans21 = ( par.transport * w->Transporters2(1) * w->C2()->Chemical(0) / (par.ka + w->C2()->Chemical(0)) );
dchem_c1[0] += trans21 - trans12;
dchem_c2[0] += trans12 - trans21;
}
};
class PIN1Localization : public WallReaction {
public:
virtual void operator()(Wall *w, double *dw1, double *dw2) {
// Cells polarize available PIN1 to Shoot Apical Meristem
if (w->C2()->BoundaryPolP()) {
if (w->AuxinSink()) {
dw1[0] = 0.; dw2[0] = 0.;
dw1[2] = 0.; dw2[2] = 0.;
// assume high auxin concentration in SAM, to convince PIN1 to polarize to it
// exocytosis regulated0
double nb_auxin = par.sam_auxin;
double receptor_level = nb_auxin * par.r / (par.kr + nb_auxin);
dw1[1] = par.k1 * w->C1()->Chemical(1) * receptor_level /( par.km + w->C1()->Chemical(1) ) - par.k2 * w->Transporters1(1);
dw2[1] = 0.;
return;
} else {
dw1[0]=dw2[0]=dw1[1]=dw2[1]=dw1[2]=dw2[2];
return;
}
}
if (w->C1()->BoundaryPolP()) {
if (w->AuxinSink()) {
dw1[0] = 0.; dw2[0] = 0.;
dw1[2] = 0.; dw2[2] = 0.;
// assume high auxin concentration in SAM, to convince PIN1 to polarize to it
// exocytosis regulated
double nb_auxin = par.sam_auxin;
double receptor_level = nb_auxin * par.r / (par.kr + nb_auxin);
dw2[1] = par.k1 * w->C2()->Chemical(1) * receptor_level /( par.km + w->C2()->Chemical(1) ) - par.k2 * w->Transporters2(1);
dw1[1] = 0.;
return;
} else {
dw1[0]=dw2[0]=dw1[1]=dw2[1]=dw1[2]=dw2[2];
return;
}
}
// PIN1 localization at wall 1
// Note: chemical 0 is Auxin (intracellular storage only)
// Chemical 1 is PIN1 (walls and intracellular storage)
//! \f$ \frac{d Pij/dt}{dt} = k_1 A_j \frac{P_i}{L_ij} - k_2 P_{ij} \f$
// Note that Pij is measured in term of concentration (mol/L)
// Pi in terms of quantity (mol)
double dPijdt1=0., dPijdt2=0.;
// normal cell
double auxin2 = w->C2()->Chemical(0);
double receptor_level1 = auxin2 * par.r / (par.kr + auxin2);
dPijdt1 =
// exocytosis regulated
par.k1 * w->C1()->Chemical(1) * receptor_level1 / ( par.km + w->C1()->Chemical(1) ) - par.k2 * w->Transporters1(1);
double auxin1 = w->C1()->Chemical(0);
double receptor_level2 = auxin1 * par.r / (par.kr + auxin1);
// normal cell
dPijdt2 =
// exocytosis regulated
par.k1 * w->C2()->Chemical(1) * receptor_level2 / ( par.km + w->C2()->Chemical(1) ) - par.k2 * w->Transporters2(1);
/* PIN1 of neighboring vascular cell inhibits PIN1 endocytosis */
dw1[0] = 0.; dw2[0] = 0.;
dw1[2] = 0.; dw2[2] = 0.;
dw1[1] = dPijdt1;
dw2[1] = dPijdt2;
}
};
inline double YlPilPli(Cell &here, Cell &nb, Wall &w) {
return nb.Chemical(2) * w.Transporters1(1) * w.Transporters2(1);
}
inline double AlPil(Cell &here, Cell &nb, Wall &w) {
return nb.Chemical(0) * w.getTransporter( &here, 1 );
}
inline double AlplusYlLil(Cell &here, Cell &nb, Wall &w) {
return ( nb.Chemical(0) + nb.Chemical(2) ) * w.Length();
}
inline double AlLil(Cell &here, Cell &nb, Wall &w) {
return ( nb.Chemical(0) ) * w.Length();
}
inline double Lil(Cell &here, Cell &nb, Wall &w) {
return w.Length();
}
inline double complex_PijAj(Cell &here, Cell &nb, Wall &w) {
// gives the amount of complex "auxinreceptor-Pin1" at the wall (at QSS)
//return here.Chemical(1) * nb.Chemical(0) / ( par.km + here.Chemical(1));
double nb_aux = (nb.BoundaryPolP() && w.AuxinSink()) ? par.sam_auxin : nb.Chemical(0);
double receptor_level = nb_aux * par.r / (par.kr + nb_aux);
return here.Chemical(1) * receptor_level / ( par.km + here.Chemical(1));
}
class AuxinAndDifferentiation : public CellReaction {
// Note: Pi and Pij measured in numbers of molecules, not concentrations
public:
virtual void operator()(Cell *c, double *dchem) {
double dPidt = 0.;
double sum_Pij = c->SumTransporters( 1 );
// exocytosis regulated:
dPidt = -par.k1 * c->ReduceCellAndWalls<double>( complex_PijAj ) + par.k2 * sum_Pij;
// production of PIN depends on auxin concentration
dPidt += (c->AtBoundaryP()?par.pin_prod_in_epidermis:par.pin_prod) * c->Chemical(0) - c->Chemical(1) * par.pin_breakdown;
/*if (c->AtBoundaryP()) {
dchem[2] = 0.01;
//cerr << "Making cell blue.\n";
} else {
dchem[2] = -0.1 * c->Chemical(2);
}*/
// no PIN production in SAM
if (c->Boundary() == Cell::SAM) {
dchem[1]=0.;
dchem[0]= - par.sam_auxin_breakdown * c->Chemical(0);
} else {
dchem[1] = dPidt;
// source of auxin
dchem[0] = par.aux_cons - par.aux_breakdown * c->Chemical(0);
}
}
};
class CellHouseKeeping {
public:
void operator() (Cell &c) const {
if (c.Boundary()==Cell::None) {
c.CheckForDivision();
// expand if this is not a provascular cell
// if (c.Chemical(2) < 0.7 ) {
c.EnlargeTargetArea((c.Chemical(0)/(1.+c.Chemical(0)))*par.cell_expansion_rate);
// }
}
}
};
void Cell::SetColor(QColor &color) {
// Red: PIN1
// Green: Auxin
// Blue: AUX1
color.setRgb(chem[1]/(1+chem[1]) * 255.,(chem[0]/(1+chem[0]) * 255.),(chem[2]/(1+chem[2]) *255.) );
}
void Cell::CheckForDivision(void) {
// if (/* Chemical(0)<0.4 && */ /* differentiated cells do not divide */ area > 2*base_area /* || Length()>50 */) {
// if (Chemical(0) > par.morphogen_div_threshold )
if (area > par.rel_cell_div_threshold * base_area ) {
/* remark no longer valid? //m->IncreaseCellCapacityIfNecessary();
// Note that calling Divide as follows prevents trouble if cell
// vector is relocated */
Divide();
}
}
// Set sum of PIN1 back to Pi_tot. Hence we assume the cell maintain a constant level of Pi.
void Cell::OnDivide(ParentInfo &parent_info, Cell &daughter) {
//cerr << "Calling Cell::OnDivide()" << endl;
// Auxin distributes between parent and daughter according to area
double area = Area(), daughter_area = daughter.Area();
double tot_area = area + daughter_area;
chem[0]*=(area/tot_area);
daughter.chem[0]*=(daughter_area/tot_area);
// For lack of detailed data, or a better rule, we assume that cells remain polarized
// after division
// So the PIN1 is redistributed according to the original polarization over the walls
// parent_info contains info about the parent
// redistribute the PIN in the endosome according to area
//chem[1] = parent_info.PINendosome*(area/tot_area);
//daughter.chem[1] = parent_info.PINendosome*(daughter_area/tot_area);
chem[1] = par.initval[1];
daughter.chem[1] = par.initval[1];
// Now redistribute the membrane PINs according to the original polarization in the parent
// mmm... I'd like to have a better, biologically motivated rule for this,
// but for lack of something better... I hope I'm excused :-). Let's say the overall
// organization of the actin fibres is not completely destroyed after division...
// distribute wallPINs according to the circumference of the parent and daughter
/* double circ = Circumference( );
double daughter_circ = daughter.Circumference();
double tot_circ = circ + daughter_circ;
double wallPINs = (circ / tot_circ) * parent_info.PINmembrane;
double daughter_wallPINs = (daughter_circ / tot_circ) * parent_info.PINmembrane;
*/
//cerr << "wallPINs = " << wallPINs << ", daughter_wallPINs = " << daughter_wallPINs << "sum = " << wallPINs + daughter_wallPINs << ", PINmembrane = " << parent_info.PINmembrane << endl;
// distrubute it according to the overall polarity
//Vector polarization = parent_info.polarization.Normalised().Perp2D();
// double sum=0.;
for (list<Wall *>::const_iterator w=walls.begin();
w!=walls.end();
w++) {
// distribute according to angle (0 degrees: maximum, 180 degrees minimum)
//double tmp=InnerProduct((*w)->getWallVector(this),polarization); // move domain from [-1,1] to [0,1]
//cerr << "[" << tmp << "]";
//sum+=tmp;
// reset transporter value
(*w)->setTransporter(this, 1, 0.);
//(*w)->setTransporter(this, 1,
}
for (list<Wall *>::const_iterator w=daughter.walls.begin();
w!=daughter.walls.end();
w++) {
// reset transporter value
(*w)->setTransporter(&daughter, 1, 0.);
}
//cerr << "Sum is " << sum << endl;
//double sum_wall_Pi = SumTransporters(1);
// After division, cells produce PIN1 (in intracellular storage) until total amount becomes Pi_tot
//SetChemical(1, par.Pi_tot - sum_wall_Pi );
//SetNewChem(1, Chemical(1));
//cerr << "[ " << sum_wall_Pi + Chemical(1) << "]";
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|