BUDEES/src/Simulateur/solver/Nucleation.h

164 lines
6.3 KiB
C

// This file is part of BUDEES, a software that allows fast simulation of bubble dynamics in electrolyzers and their interactions with electrical parameters.
//
// Copyright (C) 2024 <CNRS-ENS Rennes-Université de Rennes-CY Cergy Paris Université-Université Paris Saclay-CNAM Paris-ENS Paris-Saclay>
//
// This program is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// This program 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. See the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU General Public License along with this program. If not, see https://www.gnu.org/licenses/
/* ------------ Boundaries ------------ */
/* ------------------------------------ */
void updateNucleationSite(Parameters& parameters, MatrixXd& Sources, vector<Site>& Sites) {
int* N = parameters.N;
// Create a clone of sourcesValue since it will be modified
MatrixXd SourcesCopy(N[0], N[1]);
for (int i = 0; i < N[0]; i++) for (int j = 0; j < N[1]; j++) {
SourcesCopy(i,j) = (double) Sources(i,j);
}
// Setup delta step from hole density and check that it's not too small.
double delta = 1.0/dSite;
int di2 = round((delta / (2.0*parameters.dL[0])) - 0.5), dj2 = round((delta / (2.0*parameters.dL[1])) - 0.5); // -0.5 car conversion longueur case
di2 = MAX(1, di2); dj2 = MAX(1, dj2);
// Update all sites
for (vector<Site>::iterator it = Sites.begin(); it != Sites.end(); ++it) {
Point pos = {it->x, it->y};
int i, j; positionToIndex(N, &i, &j, pos.x/parameters.L[0], pos.y/parameters.L[1]);
double I = 0.0; int count = 0;
// Moyennage sur un rayon de 1/2d
for (int ii = -di2; ii <= +di2; ii++) for (int jj = -dj2; jj <= +dj2; jj++) {
if ((i+ii > -1) && (i+ii < N[0]) && (j+jj > -1) && (j+jj < N[1])) {
double dii = 0.0, djj = 0.0;
if (di2 > 0) dii = ((double) ii) / di2;
if (dj2 > 0) djj = ((double) jj) / dj2;
if (mathSQ(dii) + mathSQ(djj) <= 1.0) {
if (mathSign(SourcesCopy(i+ii, j+jj)) == (it->type == 0 ? +1 : -1)) {
I += ABS(SourcesCopy(i+ii, j+jj));
SourcesCopy(i+ii, j+jj) = 0;
count++;
}
}
}
}
if (I > 0.0) {
double T = Tprefactor[it->type] / I;
it->T = T;
it->growth = RFritz[it->type]/(T*growingKineticFraction);
} else it->T = INF;
}
}
vector<Site> buildNucleationSite(Parameters& parameters, Topology& topology, MatrixXd& Sources) {
int* N = parameters.N;
mathRandInit();
// Create a clone of sourcesValue since it will be modified
MatrixXd SourcesCopy(N[0], N[1]);
for (int i = 0; i < N[0]; i++) for (int j = 0; j < N[1]; j++) {
SourcesCopy(i,j) = (double) Sources(i,j);
}
// Setup delta step from hole density and check that it's not too small.
double delta = 1.0/dSite;
int di2 = round((delta / (2.0*parameters.dL[0])) - 0.5), dj2 = round((delta / (2.0*parameters.dL[1])) - 0.5); // -0.5 car conversion longueur case
if ((di2 == 0) || (dj2 == 0)) {
errorProgress("mesh refinement too small for nucleation.");
exit(1);
}
di2 = MAX(1, di2); dj2 = MAX(1, dj2);
list<Site> SiteList;
// For each poly..
for (int s = 0; s < 2; s++) for (int p = 0; p < topology.polygonCount[s]; p++) {
int n = topology.polygonIndex[s][p+1] - topology.polygonIndex[s][p];
double pos = 0.5*delta;
// For each segment...
for (int seg = 0; seg < n; seg++) {
int index = topology.polygonIndex[s][p] + seg;
int next = topology.polygonIndex[s][p] + (seg+1)%n;
Point x0 = {topology.polygonDefinitions[index].x, topology.polygonDefinitions[index].y};
Point x1 = {topology.polygonDefinitions[next].x, topology.polygonDefinitions[next].y};
// Suppression des frontière
if ( ((x0.x == 0.0) && (x1.x == 0.0)) || ((x0.x == 1.0) && (x1.x == 1.0))
|| ((x0.y == 0.0) && (x1.y == 0.0)) || ((x0.y == 1.0) && (x1.y == 1.0)) )
continue;
// DeNormalization
x0.x *= parameters.L[0]; x1.x *= parameters.L[0];
x0.y *= parameters.L[1]; x1.y *= parameters.L[1];
Vect u = {x1.x - x0.x, x1.y - x0.y};
double length = vectLength(u);
u.x = u.x / length; u.y = u.y / length;
Vect n = {u.y, -u.x};
// Parcours du segmenet
while (pos <= length) {
Point current = vectForward(x0, u, pos);
int i, j; positionToIndex(N, &i, &j, current.x/parameters.L[0], current.y/parameters.L[1]);
double I = 0.0; int count = 0;
// Moyennage sur un rayon de 1/2d
for (int ii = -di2; ii <= +di2; ii++) for (int jj = -dj2; jj <= +dj2; jj++) {
if ((i+ii > -1) && (i+ii < N[0]) && (j+jj > -1) && (j+jj < N[1])) {
double dii = 0.0, djj = 0.0;
if (di2 > 0) dii = ((double) ii) / di2;
if (dj2 > 0) djj = ((double) jj) / dj2;
if (mathSQ(dii) + mathSQ(djj) <= 1.0) {
if (mathSign(SourcesCopy(i+ii, j+jj)) == (s == 0 ? +1 : -1)) {
I += ABS(SourcesCopy(i+ii, j+jj));
SourcesCopy(i+ii, j+jj) = 0;
count++;
}
}
}
}
if (I > 0.0) {
double T = Tprefactor[s] / I;
double randomT = T*growingKineticFractionWait*mathRand();
SiteList.push_back({current.x, current.y, T, RFritz[s]/(T*growingKineticFraction), n.x, n.y, s, randomT, 0.0, 0});
}
pos += delta;
}
pos -= length;
}
}
vector<Site> Sites(SiteList.begin(), SiteList.end());
return Sites;
}
/* ------------------------------------ */
/* ----------- Nucleation ------------ */
/* ----------------------------------- */
void nucleateBubble(ListBubbles& Bubbles, vector<Site>& Sites, double t) {
if (t < Bubbles.tNextNucleation) return;
Bubbles.tNextNucleation = INF;
for (vector<Site>::iterator it = Sites.begin(); it != Sites.end(); ++it) {
if (!it->occuped) if (t >= it->t) {
// Pop a new bubble !
it->occuped = true;
it->tmiss = 0;
Bubbles.attached.push_front({it->x, it->y, 0.0, 0.0, 0.0, 1, 0, true, &(*it), {}, {0.0, 0.0}, 0.0});
Bubbles.countAttached++;
}
// May be the next one ?
if (it->t < Bubbles.tNextNucleation) {
Bubbles.tNextNucleation = it->t;
}
}
}
/* ----------------------------------- */