Line |
Branch |
Exec |
Source |
1 |
|
|
// ______ ______ _ _ _____ ______ |
2 |
|
|
// | ____| ____| | (_)/ ____| | ____| |
3 |
|
|
// | |__ | |__ | | _| (___ ___| |__ |
4 |
|
|
// | __| | __| | | | |\___ \ / __| __| |
5 |
|
|
// | | | |____| |____| |____) | (__| |____ |
6 |
|
|
// |_| |______|______|_|_____/ \___|______| |
7 |
|
|
// Finite Elements for Life Sciences and Engineering |
8 |
|
|
// |
9 |
|
|
// License: LGL2.1 License |
10 |
|
|
// FELiScE default license: LICENSE in root folder |
11 |
|
|
// |
12 |
|
|
// Main authors: E.Schenone D. Lombardi |
13 |
|
|
// |
14 |
|
|
|
15 |
|
|
// System includes |
16 |
|
|
|
17 |
|
|
// External includes |
18 |
|
|
|
19 |
|
|
// Project includes |
20 |
|
|
#include "Solver/eigenProblemPOD.hpp" |
21 |
|
|
#include "FiniteElement/elementMatrix.hpp" |
22 |
|
|
#include "FiniteElement/elementField.hpp" |
23 |
|
|
#include "Core/felisceParam.hpp" |
24 |
|
|
#include "Core/felisceTools.hpp" |
25 |
|
|
|
26 |
|
|
namespace felisce |
27 |
|
|
{ |
28 |
|
|
/*! Constructor |
29 |
|
|
*/ |
30 |
|
✗ |
EigenProblemPOD::EigenProblemPOD(): |
31 |
|
|
EigenProblem(), |
32 |
|
✗ |
m_fePotTransMemb(nullptr), |
33 |
|
✗ |
m_eigenValue(nullptr), |
34 |
|
✗ |
m_dimRomBasis(0), |
35 |
|
✗ |
m_beta(nullptr), |
36 |
|
✗ |
m_mu(nullptr), |
37 |
|
✗ |
m_gamma(nullptr), |
38 |
|
✗ |
m_eta(nullptr), |
39 |
|
✗ |
m_modeMean(nullptr), |
40 |
|
✗ |
m_fiber(nullptr), |
41 |
|
✗ |
m_source(nullptr), |
42 |
|
✗ |
m_solutionInitialized(false), |
43 |
|
✗ |
m_tensorCs(nullptr), |
44 |
|
✗ |
m_tensorE(nullptr), |
45 |
|
✗ |
m_tensorT(nullptr), |
46 |
|
✗ |
m_tensorTs(nullptr), |
47 |
|
✗ |
m_tensorY(nullptr) |
48 |
|
|
{} |
49 |
|
|
|
50 |
|
✗ |
EigenProblemPOD::~EigenProblemPOD() { |
51 |
|
✗ |
delete [] m_fePotTransMemb; |
52 |
|
✗ |
delete [] m_eigenValue; |
53 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
54 |
|
✗ |
m_basis[i].destroy(); |
55 |
|
|
} |
56 |
|
|
|
57 |
|
✗ |
m_U_0.destroy(); |
58 |
|
✗ |
m_U_0_seq.destroy(); |
59 |
|
✗ |
m_W_0.destroy(); |
60 |
|
✗ |
if (m_solutionInitialized) { |
61 |
|
✗ |
delete [] m_beta; |
62 |
|
✗ |
delete [] m_mu; |
63 |
|
✗ |
delete [] m_gamma; |
64 |
|
✗ |
delete [] m_eta; |
65 |
|
✗ |
delete [] m_modeMean; |
66 |
|
|
} |
67 |
|
✗ |
if (m_fiber) { |
68 |
|
✗ |
delete [] m_fiber; |
69 |
|
|
} |
70 |
|
✗ |
m_autocorr.destroy(); |
71 |
|
|
} |
72 |
|
|
|
73 |
|
✗ |
void EigenProblemPOD::initialize(const GeometricMeshRegion::Pointer& mesh, FelisceTransient::Pointer fstransient, MPI_Comm& comm) { |
74 |
|
✗ |
EigenProblem::initialize(mesh, fstransient, comm); |
75 |
|
|
|
76 |
|
✗ |
std::unordered_map<std::string, int> mapOfPbm; |
77 |
|
✗ |
mapOfPbm["FKPP"] = 0; |
78 |
|
✗ |
mapOfPbm["Monodomain"] = 1; |
79 |
|
✗ |
mapOfPbm["Bidomain"] = 2; |
80 |
|
✗ |
m_problem = mapOfPbm[FelisceParam::instance().model]; |
81 |
|
|
|
82 |
|
✗ |
m_dimRomBasis = FelisceParam::instance().dimRomBasis; |
83 |
|
✗ |
m_numSource = FelisceParam::instance().numberOfSource; |
84 |
|
|
|
85 |
|
|
// m_numSnapshot = FelisceParam::instance().numberOfSnapshot; |
86 |
|
|
// m_samplingFreq = FelisceParam::instance().samplingFrequency; |
87 |
|
|
|
88 |
|
✗ |
int id=0; |
89 |
|
✗ |
m_idG = id; |
90 |
|
✗ |
id++; |
91 |
|
✗ |
m_idK = id; |
92 |
|
✗ |
id++; |
93 |
|
✗ |
if (FelisceParam::instance().hasInfarct) { |
94 |
|
✗ |
m_idGs = id; |
95 |
|
✗ |
id++; |
96 |
|
|
} |
97 |
|
✗ |
m_numberOfMatrix = id; |
98 |
|
|
|
99 |
|
|
|
100 |
|
✗ |
std::vector<PhysicalVariable> listVariable; |
101 |
|
✗ |
std::vector<std::size_t> listNumComp; |
102 |
|
✗ |
if (m_problem==1) { |
103 |
|
✗ |
nameOfTheProblem = "Problem cardiac Monodomain"; |
104 |
|
✗ |
listVariable.push_back(potTransMemb); |
105 |
|
✗ |
listNumComp.push_back(1); |
106 |
|
✗ |
m_listUnknown.push_back(potTransMemb); |
107 |
|
|
} |
108 |
|
✗ |
definePhysicalVariable(listVariable,listNumComp); |
109 |
|
|
|
110 |
|
|
} |
111 |
|
|
|
112 |
|
|
// Define Physical Variable associate to the problem |
113 |
|
✗ |
void EigenProblemPOD::initPerElementType() { |
114 |
|
✗ |
if (m_problem==1) { |
115 |
|
|
//Init pointer on Finite Element, Variable or idVariable |
116 |
|
✗ |
m_ipotTransMemb = m_listVariable.getVariableIdList(potTransMemb); |
117 |
|
✗ |
m_fePotTransMemb = m_listCurrentFiniteElement[m_ipotTransMemb]; |
118 |
|
✗ |
m_elemFieldU0.initialize(DOF_FIELD,*m_fePotTransMemb); |
119 |
|
|
|
120 |
|
✗ |
if (FelisceParam::instance().hasInfarct) { |
121 |
|
✗ |
m_elemFieldFhNf0.initialize(DOF_FIELD,*m_fePotTransMemb); |
122 |
|
|
} |
123 |
|
|
|
124 |
|
|
} |
125 |
|
|
} |
126 |
|
|
|
127 |
|
|
|
128 |
|
|
|
129 |
|
|
// Assemble Matrix |
130 |
|
|
|
131 |
|
|
// mass matrix |
132 |
|
✗ |
void EigenProblemPOD::computeElementArray(const std::vector<Point*>& elemPoint, const std::vector<felInt>& elemIdPoint, felInt& iel, FlagMatrixRHS flagMatrixRHS) { |
133 |
|
|
IGNORE_UNUSED_ELEM_ID_POINT; |
134 |
|
|
IGNORE_UNUSED_FLAG_MATRIX_RHS; |
135 |
|
✗ |
if (m_problem==1) { |
136 |
|
✗ |
m_fePotTransMemb->updateMeasQuadPt(0, elemPoint); |
137 |
|
✗ |
m_fePotTransMemb->updateFirstDeriv(0, elemPoint); |
138 |
|
|
|
139 |
|
|
// Matrix[0] = phi_i * phi_j |
140 |
|
✗ |
m_elementMat[m_idG]->phi_i_phi_j(1.,*m_fePotTransMemb,0,0,1); |
141 |
|
✗ |
m_elementMat[m_idK]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor,*m_fePotTransMemb,0,0,1); |
142 |
|
|
|
143 |
|
✗ |
if (FelisceParam::instance().hasInfarct) { |
144 |
|
|
// m_Matrix[] = s * phi_i * phi_j |
145 |
|
✗ |
m_elemFieldFhNf0.setValue(m_FhNf0, *m_fePotTransMemb, iel, m_ipotTransMemb, m_ao, m_dof); |
146 |
|
✗ |
m_elementMat[m_idGs]->a_phi_i_phi_j(1.,m_elemFieldFhNf0,*m_fePotTransMemb,0,0,1); |
147 |
|
|
} |
148 |
|
|
} |
149 |
|
|
|
150 |
|
|
} |
151 |
|
|
|
152 |
|
|
|
153 |
|
|
|
154 |
|
✗ |
void EigenProblemPOD::readSnapshot(IO& io) { |
155 |
|
✗ |
m_snapshot.resize(m_numSnapshot); |
156 |
|
|
|
157 |
|
✗ |
for(int i=0; i<m_numSnapshot; i++) { |
158 |
|
|
|
159 |
|
✗ |
double* snap = nullptr; |
160 |
|
✗ |
snap = new double[m_mesh->numPoints()]; |
161 |
|
✗ |
double timeVar = 0.; |
162 |
|
|
|
163 |
|
✗ |
m_Matrix[m_idG].getVecs(m_snapshot[i],nullPetscVector); |
164 |
|
✗ |
m_snapshot[i].setFromOptions(); |
165 |
|
|
|
166 |
|
✗ |
PetscVector tmpSeq; |
167 |
|
✗ |
tmpSeq.createSeq(PETSC_COMM_SELF,m_numDof); |
168 |
|
✗ |
tmpSeq.setFromOptions(); |
169 |
|
|
|
170 |
|
|
|
171 |
|
✗ |
timeVar = i*m_samplingFreq*FelisceParam::instance().timeStep; |
172 |
|
✗ |
io.readVariable(0, timeVar,snap, m_mesh->numPoints()); |
173 |
|
|
felInt ia; |
174 |
|
✗ |
for (felInt j=0; j< m_numDof; j++) { |
175 |
|
✗ |
ia = j; |
176 |
|
✗ |
AOApplicationToPetsc(m_ao,1,&ia); |
177 |
|
✗ |
tmpSeq.setValues(1,&ia,&snap[j],INSERT_VALUES); |
178 |
|
|
} |
179 |
|
✗ |
tmpSeq.assembly(); |
180 |
|
✗ |
fromDoubleStarToVec(snap,m_snapshot[i] , m_mesh->numPoints()); |
181 |
|
|
|
182 |
|
✗ |
tmpSeq.destroy(); |
183 |
|
✗ |
delete [] snap; |
184 |
|
|
} |
185 |
|
|
|
186 |
|
|
} |
187 |
|
|
|
188 |
|
|
|
189 |
|
|
|
190 |
|
✗ |
void EigenProblemPOD::computeAutocorr() { |
191 |
|
✗ |
m_autocorr.createAIJ(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, m_numSnapshot, m_numSnapshot,PETSC_DECIDE, nullptr, PETSC_DECIDE, nullptr); |
192 |
|
✗ |
for(int i=0; i<m_numSnapshot; i++) { |
193 |
|
✗ |
for(int j=0; j<=i; j++) { |
194 |
|
✗ |
double prod = 0.0; |
195 |
|
✗ |
scalarProduct(m_snapshot[i],m_snapshot[j],prod,true); |
196 |
|
✗ |
m_autocorr.setValues(1,&i,1,&j,&prod,INSERT_VALUES); |
197 |
|
✗ |
m_autocorr.setValues(1,&j,1,&i,&prod,INSERT_VALUES); |
198 |
|
|
} |
199 |
|
|
} |
200 |
|
|
|
201 |
|
✗ |
m_autocorr.assembly(MAT_FINAL_ASSEMBLY); |
202 |
|
|
} |
203 |
|
|
|
204 |
|
|
|
205 |
|
|
|
206 |
|
|
// Visualize modes as a variable |
207 |
|
✗ |
void EigenProblemPOD::writeMode(const int iter) { |
208 |
|
|
int rankProc; |
209 |
|
✗ |
MPI_Comm_rank(m_petscComm,&rankProc); |
210 |
|
|
|
211 |
|
✗ |
std::string fileName; |
212 |
|
✗ |
if (iter > 0) |
213 |
|
✗ |
fileName = "basis" + std::to_string(iter); |
214 |
|
|
else |
215 |
|
✗ |
fileName = "basis"; |
216 |
|
|
|
217 |
|
✗ |
double* tmpSolToSave = nullptr; |
218 |
|
✗ |
tmpSolToSave = new double[m_numDof]; |
219 |
|
|
|
220 |
|
✗ |
for (int iBasis=0; iBasis<m_dimRomBasis; iBasis++) { |
221 |
|
✗ |
fromVecToDoubleStar(tmpSolToSave, m_basisPOD[iBasis], rankProc, 1); |
222 |
|
✗ |
writeEnsightVector(tmpSolToSave, iBasis, fileName); |
223 |
|
|
} |
224 |
|
✗ |
delete [] tmpSolToSave; |
225 |
|
|
|
226 |
|
✗ |
if (rankProc == 0) |
227 |
|
✗ |
writeEnsightCase(m_dimRomBasis, 1.,fileName); |
228 |
|
|
|
229 |
|
|
} |
230 |
|
|
|
231 |
|
|
|
232 |
|
|
|
233 |
|
|
|
234 |
|
|
|
235 |
|
✗ |
void EigenProblemPOD::writeEnsightSolution(const int iter) { |
236 |
|
✗ |
PetscVector solUe; |
237 |
|
✗ |
solUe.duplicateFrom(m_U_0); |
238 |
|
|
|
239 |
|
✗ |
if (iter > 0) { |
240 |
|
✗ |
projectOnDof(m_beta,m_sol,m_dimRomBasis); |
241 |
|
✗ |
if (FelisceParam::instance().model == "Bidomain") { |
242 |
|
✗ |
std::cout << "Bidomain still to be defined with POD"<<std::endl; |
243 |
|
|
} |
244 |
|
|
} else { |
245 |
|
✗ |
m_sol.copyFrom(m_U_0); |
246 |
|
✗ |
if (FelisceParam::instance().model == "Bidomain") { |
247 |
|
✗ |
std::cout << "Bidomain still to be defined with POD"<<std::endl; |
248 |
|
|
} |
249 |
|
|
} |
250 |
|
|
|
251 |
|
|
int rankProc; |
252 |
|
✗ |
MPI_Comm_rank(m_petscComm,&rankProc); |
253 |
|
|
|
254 |
|
✗ |
double* tmpSolToSave = nullptr; |
255 |
|
✗ |
tmpSolToSave = new double[m_numDof]; |
256 |
|
|
|
257 |
|
|
|
258 |
|
✗ |
if ( (FelisceParam::instance().model == "FKPP") ) { |
259 |
|
✗ |
fromVecToDoubleStar(tmpSolToSave, m_sol, rankProc, 1); |
260 |
|
✗ |
writeEnsightVector(tmpSolToSave, iter, "sol"); |
261 |
|
|
} |
262 |
|
|
|
263 |
|
✗ |
if ( (FelisceParam::instance().model == "Monodomain") || (FelisceParam::instance().model == "Bidomain") ) { |
264 |
|
✗ |
fromVecToDoubleStar(tmpSolToSave, m_sol, rankProc, 1); |
265 |
|
✗ |
writeEnsightVector(tmpSolToSave, iter, "potTransMemb"); |
266 |
|
|
} |
267 |
|
|
|
268 |
|
✗ |
if (FelisceParam::instance().model == "Bidomain") { |
269 |
|
✗ |
fromVecToDoubleStar(tmpSolToSave, solUe, rankProc, 1); |
270 |
|
✗ |
writeEnsightVector(tmpSolToSave, iter, "potExtraCell"); |
271 |
|
|
} |
272 |
|
|
|
273 |
|
✗ |
delete [] tmpSolToSave; |
274 |
|
|
|
275 |
|
✗ |
if ( (FelisceParam::instance().model == "FKPP") ) { |
276 |
|
✗ |
if (rankProc == 0) |
277 |
|
✗ |
writeEnsightCase(iter+1, FelisceParam::instance().timeStep * FelisceParam::instance().frequencyWriteSolution,"sol", FelisceParam::instance().time); |
278 |
|
|
} |
279 |
|
|
|
280 |
|
✗ |
if ( (FelisceParam::instance().model == "Monodomain") || (FelisceParam::instance().model == "Bidomain") ) { |
281 |
|
✗ |
if (rankProc == 0) |
282 |
|
✗ |
writeEnsightCase(iter+1, FelisceParam::instance().timeStep * FelisceParam::instance().frequencyWriteSolution,"potTransMemb", FelisceParam::instance().time); |
283 |
|
|
} |
284 |
|
|
|
285 |
|
✗ |
if (FelisceParam::instance().model == "Bidomain") { |
286 |
|
✗ |
std::vector<std::string> varNames; |
287 |
|
✗ |
varNames.resize(2); |
288 |
|
✗ |
varNames[0] = "potTransMemb"; |
289 |
|
✗ |
varNames[1] = "potExtraCell"; |
290 |
|
✗ |
if (rankProc == 0) |
291 |
|
✗ |
writeEnsightCase(iter+1, FelisceParam::instance().timeStep * FelisceParam::instance().frequencyWriteSolution,varNames, FelisceParam::instance().time); |
292 |
|
|
} |
293 |
|
|
|
294 |
|
|
} |
295 |
|
|
|
296 |
|
|
|
297 |
|
|
|
298 |
|
|
|
299 |
|
|
// Solve eigenproblem to build rom basis |
300 |
|
✗ |
void EigenProblemPOD::buildSolver() { |
301 |
|
✗ |
if (FelisceParam::verbose()) PetscPrintf(PETSC_COMM_WORLD,"\nBuilding Slepc Solver.\n"); |
302 |
|
|
|
303 |
|
|
#ifdef FELISCE_WITH_SLEPC |
304 |
|
✗ |
EPSSetOperators(m_eps,m_autocorr.toPetsc(),nullptr); |
305 |
|
✗ |
EPSSetDimensions(m_eps,m_numSnapshot,PETSC_DECIDE,PETSC_DECIDE); |
306 |
|
✗ |
EPSSetProblemType(m_eps,EPS_HEP); |
307 |
|
✗ |
EPSSetType(m_eps, EPSKRYLOVSCHUR); //EPSARNOLDI); |
308 |
|
✗ |
EPSSetTolerances(m_eps,1.0e-5,200); |
309 |
|
✗ |
EPSKrylovSchurSetRestart(m_eps,0.6); |
310 |
|
✗ |
EPSSetFromOptions(m_eps); |
311 |
|
|
#endif |
312 |
|
|
} |
313 |
|
|
|
314 |
|
|
|
315 |
|
✗ |
void EigenProblemPOD::solve() { |
316 |
|
✗ |
if (FelisceParam::verbose()) PetscPrintf(PETSC_COMM_WORLD,"\nSolve eigenvalue problem.\n"); |
317 |
|
|
#ifdef FELISCE_WITH_SLEPC |
318 |
|
|
PetscInt nconv; |
319 |
|
|
//First request a singular value from one end of the spectrum |
320 |
|
✗ |
EPSSetWhichEigenpairs(m_eps,EPS_LARGEST_REAL); |
321 |
|
✗ |
if (m_verbose >= 10) { |
322 |
|
✗ |
m_autocorr.view(); |
323 |
|
✗ |
EPSView(m_eps,PETSC_VIEWER_STDOUT_WORLD); |
324 |
|
|
} |
325 |
|
|
|
326 |
|
✗ |
EPSSolve(m_eps); |
327 |
|
|
|
328 |
|
|
//Get number of converged singular values |
329 |
|
✗ |
EPSGetConverged(m_eps,&nconv); |
330 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD," Number of converged values = %d.\n", nconv); |
331 |
|
|
|
332 |
|
✗ |
if ( m_dimRomBasis > m_numDof ) { |
333 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD," Warning! Number of function basis bigger than number of dof.\n"); |
334 |
|
✗ |
m_dimRomBasis = m_numDof; |
335 |
|
|
} |
336 |
|
|
|
337 |
|
✗ |
if ( m_dimRomBasis > nconv ) { |
338 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD," Warning! Number of function basis bigger than number of converged values.\n"); |
339 |
|
✗ |
m_dimRomBasis = nconv; |
340 |
|
|
} |
341 |
|
|
|
342 |
|
|
double sigma; |
343 |
|
|
//Get converged singular values: largest singular value is stored in sigma. |
344 |
|
✗ |
if (nconv > 0) { |
345 |
|
✗ |
for (int i=0; i < m_dimRomBasis; i++) { |
346 |
|
✗ |
EPSGetEigenpair(m_eps,i,&sigma,FELISCE_PETSC_NULLPTR,FELISCE_PETSC_NULLPTR,FELISCE_PETSC_NULLPTR); |
347 |
|
|
} |
348 |
|
|
} else { |
349 |
|
✗ |
m_dimRomBasis = 0; |
350 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD," Unable to compute smallest singular value!\n\n"); |
351 |
|
✗ |
exit(1); |
352 |
|
|
} |
353 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "m_dimRomBasis = %d \n", m_dimRomBasis); |
354 |
|
|
|
355 |
|
✗ |
if (m_eigenValue == nullptr) |
356 |
|
✗ |
m_eigenValue = new double[m_dimRomBasis]; |
357 |
|
|
|
358 |
|
|
//m_basis.clear(); |
359 |
|
|
//m_basis.resize(m_dimRomBasis); |
360 |
|
✗ |
m_basisPOD.resize(m_dimRomBasis); |
361 |
|
✗ |
m_eigenvector.resize(m_dimRomBasis); |
362 |
|
|
double norm2; |
363 |
|
|
|
364 |
|
✗ |
double imgpart = 0.; |
365 |
|
✗ |
for (int i=0; i < m_dimRomBasis; i++) { |
366 |
|
✗ |
m_autocorr.getVecs(nullPetscVector,m_eigenvector[i]); |
367 |
|
✗ |
EPSGetEigenpair(m_eps,i,&m_eigenValue[i],&imgpart,m_eigenvector[i].toPetsc(),FELISCE_PETSC_NULLPTR); |
368 |
|
✗ |
if (Tools::notEqual(imgpart,0.0)) { |
369 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "Warning! Im(m_eigenValue[%d]) = %e \n", i, imgpart); |
370 |
|
|
} |
371 |
|
✗ |
m_eigenvector[i].norm(NORM_2,&norm2); |
372 |
|
|
|
373 |
|
✗ |
if (m_verbose >= 1) { |
374 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "m_eigenValue[%d] = %e \n", i, m_eigenValue[i]); |
375 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "norm2[%d] = %e \n", i, norm2); |
376 |
|
|
} |
377 |
|
|
} |
378 |
|
|
|
379 |
|
|
felInt ia; |
380 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
381 |
|
✗ |
m_Matrix[0].getVecs(nullPetscVector,m_basisPOD[i]); |
382 |
|
✗ |
for(int h=0; h<m_numSnapshot; h++) { |
383 |
|
✗ |
ia = h; |
384 |
|
✗ |
AOApplicationToPetsc(m_ao,1,&ia); |
385 |
|
✗ |
double tmpCoeff = 0.0; |
386 |
|
✗ |
m_eigenvector[i].getValues(1,&ia,&tmpCoeff); |
387 |
|
✗ |
m_basisPOD[i].axpy( tmpCoeff, m_snapshot[h]); |
388 |
|
|
}; |
389 |
|
✗ |
double tmpRescale = m_eigenValue[i]; |
390 |
|
✗ |
tmpRescale = std::sqrt(tmpRescale); |
391 |
|
✗ |
tmpRescale = 1.0/tmpRescale; |
392 |
|
✗ |
m_basisPOD[i].scale(tmpRescale); |
393 |
|
|
}; |
394 |
|
|
|
395 |
|
✗ |
std::string fileName = FelisceParam::instance().resultDir + "/eigenValue"; |
396 |
|
✗ |
viewPOD(m_eigenValue,m_dimRomBasis,fileName); |
397 |
|
|
#endif |
398 |
|
|
} |
399 |
|
|
|
400 |
|
|
|
401 |
|
|
|
402 |
|
|
|
403 |
|
|
// Reduced model functions |
404 |
|
|
//======================== |
405 |
|
|
// Initialization |
406 |
|
✗ |
void EigenProblemPOD::readData(IO& io) { |
407 |
|
✗ |
m_Matrix[0].getVecs(m_U_0,nullPetscVector); |
408 |
|
✗ |
m_U_0.setFromOptions(); |
409 |
|
|
|
410 |
|
✗ |
m_W_0.duplicateFrom(m_U_0); |
411 |
|
✗ |
m_W_0.copyFrom(m_U_0); |
412 |
|
|
|
413 |
|
✗ |
m_U_0_seq.createSeq(PETSC_COMM_SELF,m_numDof); |
414 |
|
✗ |
m_U_0_seq.setFromOptions(); |
415 |
|
✗ |
if (FelisceParam::instance().hasInitialCondition) { |
416 |
|
|
// Read initial solution file (*.00000.scl) |
417 |
|
✗ |
double* initialSolution = nullptr; |
418 |
|
✗ |
initialSolution = new double[m_mesh->numPoints()]; |
419 |
|
✗ |
double* ionicSolution = nullptr; |
420 |
|
✗ |
ionicSolution = new double[m_mesh->numPoints()]; |
421 |
|
|
// potTransMemb |
422 |
|
✗ |
io.readVariable(0, 0.,initialSolution, m_mesh->numPoints()); |
423 |
|
|
// ionicVar |
424 |
|
✗ |
io.readVariable(1, 0.,ionicSolution, m_mesh->numPoints()); |
425 |
|
|
felInt ia; |
426 |
|
✗ |
for (felInt i=0; i< m_numDof; i++) { |
427 |
|
✗ |
ia = i; |
428 |
|
✗ |
AOApplicationToPetsc(m_ao,1,&ia); |
429 |
|
✗ |
m_U_0_seq.setValues(1,&ia,&initialSolution[i],INSERT_VALUES); |
430 |
|
|
} |
431 |
|
✗ |
m_U_0_seq.assembly(); |
432 |
|
✗ |
fromDoubleStarToVec(initialSolution, m_U_0, m_mesh->numPoints()); |
433 |
|
✗ |
fromDoubleStarToVec(ionicSolution, m_W_0, m_mesh->numPoints()); |
434 |
|
✗ |
delete [] initialSolution; |
435 |
|
✗ |
delete [] ionicSolution; |
436 |
|
|
// Read fibers file (*.00000.vct) |
437 |
|
✗ |
if (FelisceParam::instance().testCase == 1) { |
438 |
|
✗ |
if (m_fiber == nullptr) |
439 |
|
✗ |
m_fiber = new double[m_mesh->numPoints()*3]; |
440 |
|
✗ |
io.readVariable(2, 0.,m_fiber, m_mesh->numPoints()*3); |
441 |
|
|
} |
442 |
|
|
} else { |
443 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "U_0 not read from file (zero vector)."); |
444 |
|
✗ |
m_U_0.set( 0.); |
445 |
|
✗ |
m_U_0.assembly(); |
446 |
|
✗ |
m_W_0.set( 0.); |
447 |
|
✗ |
m_W_0.assembly(); |
448 |
|
✗ |
m_U_0_seq.set( 0.); |
449 |
|
✗ |
m_U_0_seq.assembly(); |
450 |
|
|
} |
451 |
|
|
} |
452 |
|
|
|
453 |
|
✗ |
void EigenProblemPOD::getFiberDirection(felInt iel, int iUnknown, std::vector<double>& elemFiber) { |
454 |
|
✗ |
int numSupport = static_cast<int>(m_supportDofUnknown[iUnknown].iEle()[iel+1] - m_supportDofUnknown[iUnknown].iEle()[iel]); |
455 |
|
✗ |
elemFiber.resize( numSupport*3,0.); |
456 |
|
✗ |
for (int i = 0; i < numSupport; i++) { |
457 |
|
✗ |
elemFiber[i*3] = m_fiber[3*(m_supportDofUnknown[iUnknown].iSupportDof()[m_supportDofUnknown[iUnknown].iEle()[iel]+i])]; |
458 |
|
✗ |
elemFiber[i*3+1] = m_fiber[3*(m_supportDofUnknown[iUnknown].iSupportDof()[m_supportDofUnknown[iUnknown].iEle()[iel]+i])+1]; |
459 |
|
✗ |
elemFiber[i*3+2] = m_fiber[3*(m_supportDofUnknown[iUnknown].iSupportDof()[m_supportDofUnknown[iUnknown].iEle()[iel]+i])+2]; |
460 |
|
|
} |
461 |
|
|
} |
462 |
|
|
|
463 |
|
✗ |
void EigenProblemPOD::initializeROM() { |
464 |
|
✗ |
if (FelisceParam::verbose() > 10 ) PetscPrintf(PETSC_COMM_WORLD,"\nEigenProblemPOD::initializeROM()\n"); |
465 |
|
|
|
466 |
|
✗ |
m_basisPOD.resize(m_dimRomBasis); |
467 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
468 |
|
✗ |
m_Matrix[m_idG].getVecs(m_basisPOD[i],nullPetscVector); |
469 |
|
✗ |
m_basisPOD[i].setFromOptions(); |
470 |
|
✗ |
double* vec = new double[m_numDof]; |
471 |
|
✗ |
readEnsightFile(vec,"basis",i,m_numDof); |
472 |
|
✗ |
fromDoubleStarToVec(vec, m_basisPOD[i], m_numDof); |
473 |
|
✗ |
delete [] vec; |
474 |
|
✗ |
if (FelisceParam::verbose() > 40) m_basisPOD[i].view(); |
475 |
|
|
} |
476 |
|
|
|
477 |
|
✗ |
readEigenValueFromFile(); |
478 |
|
|
|
479 |
|
✗ |
initializeSolution(); |
480 |
|
|
} |
481 |
|
|
|
482 |
|
✗ |
void EigenProblemPOD::readEnsightFile(double* vec, std::string varName, int idIter, felInt size) { |
483 |
|
✗ |
std::string iteration; |
484 |
|
✗ |
if (idIter < 10) |
485 |
|
✗ |
iteration = "0000" + std::to_string(idIter); |
486 |
|
✗ |
else if (idIter < 100) |
487 |
|
✗ |
iteration = "000" + std::to_string(idIter); |
488 |
|
✗ |
else if (idIter < 1000) |
489 |
|
✗ |
iteration = "00" + std::to_string(idIter); |
490 |
|
✗ |
else if (idIter < 10000) |
491 |
|
✗ |
iteration = "0" + std::to_string(idIter); |
492 |
|
✗ |
else if (idIter < 100000) |
493 |
|
✗ |
iteration = std::to_string(idIter); |
494 |
|
|
|
495 |
|
✗ |
std::string fileName = FelisceParam::instance().inputDirectory + "/" + varName + "." + iteration + ".scl"; |
496 |
|
✗ |
if (FelisceParam::verbose() > 1) PetscPrintf(PETSC_COMM_WORLD, "Read file %s\n", fileName.c_str()); |
497 |
|
✗ |
FILE* pFile = fopen (fileName.c_str(),"r"); |
498 |
|
✗ |
if (pFile==nullptr) { |
499 |
|
✗ |
FEL_ERROR("Impossible to read "+fileName +"."); |
500 |
|
|
} |
501 |
|
|
char str[80]; |
502 |
|
✗ |
std::string display; |
503 |
|
✗ |
if ( fscanf(pFile,"%s", str) !=1 ) { |
504 |
|
✗ |
FEL_WARNING("Failed ot read str"); |
505 |
|
|
} |
506 |
|
✗ |
display = str; |
507 |
|
✗ |
if ( fscanf(pFile,"%s", str) !=1 ) { |
508 |
|
✗ |
FEL_WARNING("Failed ot read str"); |
509 |
|
|
} |
510 |
|
✗ |
display = str; |
511 |
|
✗ |
if ( fscanf(pFile,"%s", str) !=1) { |
512 |
|
✗ |
FEL_WARNING("Failed ot read str"); |
513 |
|
|
} |
514 |
|
✗ |
display = str; |
515 |
|
|
//Read values |
516 |
|
✗ |
int count = 0; |
517 |
|
✗ |
for ( felInt i = 0; i < size; i++) { |
518 |
|
✗ |
if ( fscanf(pFile,"%lf", &vec[i]) !=1 ) { |
519 |
|
✗ |
FEL_WARNING("Failed ot read vec"); |
520 |
|
|
} |
521 |
|
✗ |
if ( count == 6 ) { |
522 |
|
✗ |
if ( fscanf( pFile, "\n" ) !=0 ) { |
523 |
|
✗ |
FEL_WARNING("Failed ot read"); |
524 |
|
|
} |
525 |
|
✗ |
count=0; |
526 |
|
|
} |
527 |
|
✗ |
count++; |
528 |
|
|
} |
529 |
|
✗ |
fclose(pFile); |
530 |
|
|
} |
531 |
|
|
|
532 |
|
✗ |
void EigenProblemPOD::readEigenValueFromFile() { |
533 |
|
✗ |
std::string fileName = FelisceParam::instance().inputDirectory + "/eigenValue"; |
534 |
|
✗ |
if (FelisceParam::verbose() > 1) PetscPrintf(PETSC_COMM_WORLD, "Read file %s\n", fileName.c_str()); |
535 |
|
✗ |
FILE* pFile = fopen (fileName.c_str(),"r"); |
536 |
|
✗ |
if (pFile==nullptr) { |
537 |
|
✗ |
FEL_ERROR("Impossible to read "+fileName +"."); |
538 |
|
|
} |
539 |
|
|
//Read values |
540 |
|
✗ |
m_eigenValue = new double[m_dimRomBasis]; |
541 |
|
✗ |
for ( felInt i = 0; i < m_dimRomBasis; i++) { |
542 |
|
✗ |
if ( fscanf(pFile,"%lf", &m_eigenValue[i]) !=1 ) { |
543 |
|
✗ |
FEL_WARNING("Failed ot read vec"); |
544 |
|
|
} |
545 |
|
✗ |
if ( fscanf( pFile, "\n" ) !=0 ) { |
546 |
|
✗ |
FEL_WARNING("Failed ot read"); |
547 |
|
|
} |
548 |
|
|
} |
549 |
|
✗ |
fclose(pFile); |
550 |
|
|
|
551 |
|
✗ |
if (FelisceParam::verbose() > 10) { |
552 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "EigenValues:\n"); |
553 |
|
✗ |
for ( felInt i = 0; i < m_dimRomBasis; i++) { |
554 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "[%d] %e \n", i, m_eigenValue[i]); |
555 |
|
|
} |
556 |
|
|
} |
557 |
|
|
|
558 |
|
|
} |
559 |
|
|
|
560 |
|
|
|
561 |
|
✗ |
void EigenProblemPOD::computeGamma() { |
562 |
|
✗ |
if (FelisceParam::verbose() > 10 ) PetscPrintf(PETSC_COMM_WORLD,"\nEigenProblemPOD::computeGamma()\n"); |
563 |
|
|
|
564 |
|
|
// Warning: for Monodomain heart potential problem coupled with FhN |
565 |
|
|
// todo To be generalized |
566 |
|
|
|
567 |
|
✗ |
double Am = FelisceParam::instance().Am; |
568 |
|
✗ |
double Cm = FelisceParam::instance().Cm; |
569 |
|
✗ |
double s = FelisceParam::instance().f0; |
570 |
|
✗ |
double a = FelisceParam::instance().alpha; |
571 |
|
✗ |
double eps = FelisceParam::instance().epsilon; |
572 |
|
✗ |
double gam = FelisceParam::instance().gammaEl; |
573 |
|
|
|
574 |
|
✗ |
if (m_gamma == nullptr) |
575 |
|
✗ |
m_gamma = new double[m_dimRomBasis]; |
576 |
|
✗ |
if (m_eta == nullptr) |
577 |
|
✗ |
m_eta = new double[m_dimRomBasis]; |
578 |
|
|
|
579 |
|
✗ |
for (int p=0; p<m_dimRomBasis; p++) { |
580 |
|
✗ |
m_gamma[p] = 0.0; |
581 |
|
|
}; |
582 |
|
|
|
583 |
|
✗ |
if ( !FelisceParam::instance().hasInfarct ) { |
584 |
|
✗ |
for (int p=0; p<m_dimRomBasis; p++) { |
585 |
|
✗ |
m_gamma[p] = - ( s*a*Am ) * m_beta[p]; |
586 |
|
✗ |
m_gamma[p] += - Am * m_mu[p]; |
587 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
588 |
|
✗ |
m_gamma[p] -= m_beta[i]*m_tensorE[secondOrderGlobalIndex(p,i)]; |
589 |
|
✗ |
for (int j=0; j<m_dimRomBasis; j++) { |
590 |
|
✗ |
m_gamma[p] += ( s*(a+1)*Am) * m_beta[i] * m_beta[j] * m_tensorT[thirdOrderGlobalIndex(i,j,p)]; |
591 |
|
✗ |
for (int k=0; k<m_dimRomBasis; k++) { |
592 |
|
✗ |
m_gamma[p] += - s*Am * m_beta[i] * m_beta[j] * m_beta[k] * m_tensorY[fourthOrderGlobalIndex(i,j,k,p)]; |
593 |
|
|
} |
594 |
|
|
} |
595 |
|
|
} |
596 |
|
|
|
597 |
|
✗ |
if (FelisceParam::instance().hasSource) { |
598 |
|
✗ |
double delay[m_numSource]; |
599 |
|
✗ |
for (int iS=0; iS<m_numSource; iS++) { |
600 |
|
✗ |
delay[iS] = FelisceParam::instance().delayStim[iS]; |
601 |
|
|
} |
602 |
|
✗ |
double tPeriod=fmod(m_fstransient->time,FelisceParam::instance().timePeriod); |
603 |
|
✗ |
for (int iS=0; iS<m_numSource; iS++) { |
604 |
|
✗ |
if ((tPeriod>=delay[iS])&&(tPeriod<=(delay[iS]+FelisceParam::instance().stimTime))) |
605 |
|
✗ |
m_gamma[p] += Am*m_source[iS][p]; |
606 |
|
|
} |
607 |
|
|
} |
608 |
|
|
|
609 |
|
✗ |
m_gamma[p] = m_gamma[p] / (Am*Cm); |
610 |
|
|
|
611 |
|
✗ |
m_eta[p] = eps * ( gam * m_beta[p] - m_mu[p] ); |
612 |
|
|
|
613 |
|
|
} |
614 |
|
✗ |
} else if (FelisceParam::instance().hasInfarct) { |
615 |
|
✗ |
for (int p=0; p<m_dimRomBasis; p++) { |
616 |
|
✗ |
m_gamma[p] += - Am * m_mu[p]; |
617 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
618 |
|
✗ |
m_gamma[p] -= m_beta[i]*m_tensorE[secondOrderGlobalIndex(p,i)]; |
619 |
|
✗ |
m_gamma[p] += - a*Am * m_beta[i] * m_tensorCs[secondOrderGlobalIndex(i,p)]; |
620 |
|
✗ |
for (int j=0; j<m_dimRomBasis; j++) { |
621 |
|
✗ |
m_gamma[p] += (a+1)*Am * m_beta[i] * m_beta[j] * m_tensorTs[thirdOrderGlobalIndex(i,j,p)]; |
622 |
|
✗ |
for (int k=0; k<m_dimRomBasis; k++) { |
623 |
|
✗ |
m_gamma[p] += - Am * m_beta[i] * m_beta[j] * m_beta[k] * m_tensorY[fourthOrderGlobalIndex(i,j,k,p)]; |
624 |
|
|
} |
625 |
|
|
} |
626 |
|
|
} |
627 |
|
|
|
628 |
|
✗ |
if (FelisceParam::instance().hasSource) { |
629 |
|
✗ |
double delay[m_numSource]; |
630 |
|
✗ |
for (int iS=0; iS<m_numSource; iS++) { |
631 |
|
✗ |
delay[iS] = FelisceParam::instance().delayStim[iS]; |
632 |
|
|
} |
633 |
|
✗ |
double tPeriod=fmod(m_fstransient->time,FelisceParam::instance().timePeriod); |
634 |
|
✗ |
for (int iS=0; iS<m_numSource; iS++) { |
635 |
|
✗ |
if ((tPeriod>=delay[iS])&&(tPeriod<=(delay[iS]+FelisceParam::instance().stimTime))) |
636 |
|
✗ |
m_gamma[p] += Am*m_source[iS][p]; |
637 |
|
|
} |
638 |
|
|
} |
639 |
|
|
|
640 |
|
✗ |
m_gamma[p] = m_gamma[p] / (Am*Cm); |
641 |
|
|
|
642 |
|
✗ |
m_eta[p] = eps * ( gam * m_beta[p] - m_mu[p] ); |
643 |
|
|
|
644 |
|
|
} |
645 |
|
|
} |
646 |
|
|
|
647 |
|
|
} |
648 |
|
|
|
649 |
|
|
|
650 |
|
|
|
651 |
|
|
|
652 |
|
|
|
653 |
|
✗ |
void EigenProblemPOD::initializeSolution() { |
654 |
|
✗ |
if (FelisceParam::verbose() > 10 ) PetscPrintf(PETSC_COMM_WORLD,"\nEigenProblemPOD::initializeSolution()\n"); |
655 |
|
|
|
656 |
|
✗ |
if ( !m_solutionInitialized) { |
657 |
|
✗ |
if (m_beta == nullptr) |
658 |
|
✗ |
m_beta = new double[m_dimRomBasis]; |
659 |
|
|
|
660 |
|
✗ |
projectOnBasis(m_U_0,m_beta,true); |
661 |
|
|
|
662 |
|
✗ |
if (m_mu == nullptr) |
663 |
|
✗ |
m_mu = new double[m_dimRomBasis]; |
664 |
|
|
|
665 |
|
✗ |
projectOnBasis(m_W_0,m_mu,true); |
666 |
|
|
|
667 |
|
✗ |
m_solutionInitialized = true; |
668 |
|
|
} |
669 |
|
|
} |
670 |
|
|
|
671 |
|
|
|
672 |
|
✗ |
void EigenProblemPOD::setFhNf0(std::vector<double>& valuef0) { |
673 |
|
✗ |
m_FhNf0.createSeq(PETSC_COMM_SELF,m_numDof); |
674 |
|
✗ |
m_FhNf0.setFromOptions(); |
675 |
|
|
|
676 |
|
|
// initialize SEQUENTIAL vectors (in order to use elemField) |
677 |
|
|
felInt ia; |
678 |
|
✗ |
for (felInt i=0; i< m_numDof; i++) { |
679 |
|
✗ |
ia = i; |
680 |
|
✗ |
AOApplicationToPetsc(m_ao,1,&ia); |
681 |
|
✗ |
m_FhNf0.setValues(1,&ia,&valuef0[i],INSERT_VALUES); |
682 |
|
|
} |
683 |
|
✗ |
m_FhNf0.assembly(); |
684 |
|
|
|
685 |
|
|
} |
686 |
|
|
|
687 |
|
✗ |
void EigenProblemPOD::setIapp(std::vector<double>& iApp, int& idS) { |
688 |
|
✗ |
if (m_source == nullptr) { |
689 |
|
✗ |
m_source = new double*[m_numSource]; |
690 |
|
✗ |
for (int i=0; i<m_numSource; i++) { |
691 |
|
✗ |
m_source[i] = new double[m_dimRomBasis]; |
692 |
|
|
} |
693 |
|
|
} |
694 |
|
|
|
695 |
|
✗ |
PetscVector RHS; |
696 |
|
✗ |
RHS.createSeq(PETSC_COMM_SELF,m_numDof); |
697 |
|
✗ |
RHS.setFromOptions(); |
698 |
|
|
// initialize SEQUENTIAL vectors (in order to use elemField) |
699 |
|
|
felInt ia; |
700 |
|
✗ |
for (felInt i=0; i< m_numDof; i++) { |
701 |
|
✗ |
ia = i; |
702 |
|
✗ |
AOApplicationToPetsc(m_ao,1,&ia); |
703 |
|
✗ |
RHS.setValues(1,&ia,&iApp[i],INSERT_VALUES); |
704 |
|
|
} |
705 |
|
✗ |
RHS.assembly(); |
706 |
|
|
|
707 |
|
✗ |
projectOnBasis(RHS, m_source[idS], true); |
708 |
|
|
|
709 |
|
✗ |
RHS.destroy(); |
710 |
|
|
} |
711 |
|
|
|
712 |
|
|
|
713 |
|
|
|
714 |
|
|
|
715 |
|
|
|
716 |
|
|
|
717 |
|
✗ |
void EigenProblemPOD::checkSolution(PetscVector& solution, PetscVector& refSol, double & error) { |
718 |
|
|
|
719 |
|
✗ |
PetscVector tmpVec1; |
720 |
|
✗ |
tmpVec1.duplicateFrom(refSol); |
721 |
|
|
|
722 |
|
✗ |
PetscVector tmpVec2; |
723 |
|
✗ |
tmpVec2.duplicateFrom(refSol); |
724 |
|
|
|
725 |
|
✗ |
waxpy(tmpVec1,-1.0,refSol,solution); |
726 |
|
|
|
727 |
|
✗ |
mult(m_Matrix[m_idG],tmpVec1,tmpVec2); |
728 |
|
|
|
729 |
|
✗ |
dot(tmpVec1,tmpVec2,&error); |
730 |
|
|
|
731 |
|
✗ |
tmpVec1.destroy(); |
732 |
|
✗ |
tmpVec2.destroy(); |
733 |
|
|
} |
734 |
|
|
|
735 |
|
|
|
736 |
|
|
|
737 |
|
✗ |
void EigenProblemPOD::checkBasis() { |
738 |
|
|
|
739 |
|
✗ |
if (!m_solutionInitialized) { |
740 |
|
✗ |
initializeSolution(); |
741 |
|
|
} |
742 |
|
|
|
743 |
|
✗ |
double convergence[m_dimRomBasis]; |
744 |
|
|
int rankProc; |
745 |
|
✗ |
MPI_Comm_rank(m_petscComm,&rankProc); |
746 |
|
✗ |
PetscVector tmpSol; |
747 |
|
✗ |
tmpSol.duplicateFrom(m_snapshot[m_numSnapshot-1]); |
748 |
|
|
|
749 |
|
|
// Check convergence on the last snapshot |
750 |
|
✗ |
double* tmpSolToSave = nullptr; |
751 |
|
✗ |
tmpSolToSave = new double[m_numDof]; |
752 |
|
✗ |
for (int iBasis=0; iBasis < m_dimRomBasis; iBasis++) { |
753 |
|
✗ |
double tmpCoeff[iBasis]; |
754 |
|
✗ |
for (int jBasis=0; jBasis<iBasis+1; jBasis++) { |
755 |
|
✗ |
tmpCoeff[jBasis] = 0.0; |
756 |
|
✗ |
scalarProduct(m_snapshot[m_numSnapshot-1], m_basisPOD[jBasis], tmpCoeff[jBasis],true); |
757 |
|
|
} |
758 |
|
|
|
759 |
|
✗ |
projectOnDof(tmpCoeff,tmpSol,iBasis); |
760 |
|
|
|
761 |
|
✗ |
checkSolution(tmpSol,m_snapshot[m_numSnapshot-1],convergence[iBasis]); |
762 |
|
|
|
763 |
|
✗ |
fromVecToDoubleStar(tmpSolToSave, tmpSol, rankProc, 1); |
764 |
|
✗ |
writeEnsightVector(tmpSolToSave, iBasis, "potTransMembConv"); |
765 |
|
|
|
766 |
|
✗ |
PetscPrintf(PETSC_COMM_WORLD, "convergence[%d] = %e \n", iBasis, convergence[iBasis]); |
767 |
|
|
} |
768 |
|
✗ |
if (rankProc == 0) |
769 |
|
✗ |
writeEnsightCase(m_dimRomBasis, 1., "potTransMembConv"); |
770 |
|
✗ |
std::string fileName = FelisceParam::instance().resultDir + "/convergence"; |
771 |
|
✗ |
viewPOD(convergence,m_dimRomBasis,fileName); |
772 |
|
|
|
773 |
|
✗ |
delete [] tmpSolToSave; |
774 |
|
✗ |
tmpSol.destroy(); |
775 |
|
|
} |
776 |
|
|
|
777 |
|
|
|
778 |
|
✗ |
void EigenProblemPOD::checkOrthonormality() { |
779 |
|
✗ |
std::ofstream fOrth; |
780 |
|
✗ |
fOrth.open("Orthonormality.log"); |
781 |
|
|
|
782 |
|
✗ |
std::cout << "Checking basis orthonormality: " <<std::endl; |
783 |
|
✗ |
for(int i=0; i<m_dimRomBasis; i++) { |
784 |
|
✗ |
for(int j=0; j<i+1; j++) { |
785 |
|
✗ |
double prod = 0.0; |
786 |
|
✗ |
scalarProduct(m_basisPOD[i], m_basisPOD[j], prod, true); |
787 |
|
✗ |
fOrth << "Mass("<<i<<","<<j<<") = " << prod <<std::endl; |
788 |
|
|
}; |
789 |
|
|
}; |
790 |
|
✗ |
fOrth.close(); |
791 |
|
|
} |
792 |
|
|
|
793 |
|
|
|
794 |
|
|
|
795 |
|
✗ |
void EigenProblemPOD::scalarProduct(PetscVector& firstVec, PetscVector& secondVec, double& prod, bool flagMass) { |
796 |
|
✗ |
PetscVector tmpVec; |
797 |
|
✗ |
tmpVec.duplicateFrom(firstVec); |
798 |
|
|
|
799 |
|
✗ |
if (flagMass) { |
800 |
|
✗ |
mult(m_Matrix[m_idG], firstVec, tmpVec); |
801 |
|
✗ |
dot(tmpVec,secondVec,&prod); |
802 |
|
|
} else { |
803 |
|
✗ |
dot( firstVec,secondVec,&prod); |
804 |
|
|
} |
805 |
|
|
|
806 |
|
✗ |
tmpVec.destroy(); |
807 |
|
✗ |
}; |
808 |
|
|
|
809 |
|
|
|
810 |
|
|
|
811 |
|
|
// Modified Gram-Schmidt (see Golub & Van Loan) |
812 |
|
✗ |
void EigenProblemPOD::MGS(std::vector<PetscVector>& subspace) { |
813 |
|
✗ |
PetscVector tmpVec; |
814 |
|
✗ |
tmpVec.duplicateFrom(subspace[0]); |
815 |
|
✗ |
int dimSubSpace = subspace.size(); |
816 |
|
|
|
817 |
|
✗ |
for(int iBasis=0; iBasis<dimSubSpace; iBasis++) { |
818 |
|
✗ |
mult(m_Matrix[0],subspace[iBasis],tmpVec); |
819 |
|
✗ |
double rii = 0.0; |
820 |
|
✗ |
dot(tmpVec,subspace[iBasis],&rii); |
821 |
|
✗ |
rii = 1.0/std::sqrt(rii); |
822 |
|
✗ |
subspace[iBasis].scale(rii); |
823 |
|
|
|
824 |
|
✗ |
mult(m_Matrix[0],subspace[iBasis],tmpVec); |
825 |
|
|
|
826 |
|
✗ |
for(int j=iBasis+1; j<m_dimRomBasis; j++) { |
827 |
|
✗ |
double rij = 0.0; |
828 |
|
✗ |
dot(tmpVec,subspace[j],&rij); |
829 |
|
✗ |
rij = -1.0 * rij; |
830 |
|
✗ |
subspace[j].axpy(rij,subspace[iBasis]); |
831 |
|
|
} |
832 |
|
|
} |
833 |
|
✗ |
tmpVec.destroy(); |
834 |
|
|
} |
835 |
|
|
|
836 |
|
|
|
837 |
|
|
|
838 |
|
✗ |
void EigenProblemPOD::computeMatrixZeta() { |
839 |
|
✗ |
if (FelisceParam::verbose() > 10 ) PetscPrintf(PETSC_COMM_WORLD,"\nEigenProblemPOD::computeMatrixZeta()\n"); |
840 |
|
✗ |
if ( m_zeta.size() == 0 ) { |
841 |
|
✗ |
m_size2 = static_cast<int>(m_dimRomBasis*(m_dimRomBasis+1)/2); |
842 |
|
✗ |
m_zeta.resize(m_size2); |
843 |
|
✗ |
for (int i=0; i<m_size2; i++) |
844 |
|
✗ |
m_zeta[i].duplicateFrom(m_basisPOD[0]); |
845 |
|
|
} |
846 |
|
|
|
847 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { // row |
848 |
|
✗ |
for (int j=0; j<i+1; j++) { // column |
849 |
|
|
// secondOrderGlobalIndex(i,j) = i*(i+1)/2 + j |
850 |
|
✗ |
pointwiseMult(m_zeta[secondOrderGlobalIndex(i,j)],m_basisPOD[i],m_basisPOD[j]); |
851 |
|
|
} |
852 |
|
|
} |
853 |
|
|
|
854 |
|
|
} |
855 |
|
|
|
856 |
|
|
|
857 |
|
|
|
858 |
|
✗ |
int EigenProblemPOD::secondOrderGlobalIndex(int& i, int& j) { |
859 |
|
|
int id; |
860 |
|
|
|
861 |
|
✗ |
if (i >= j) { |
862 |
|
✗ |
id = static_cast<int>(i*(i+1)/2+j); |
863 |
|
|
} else { |
864 |
|
✗ |
id = static_cast<int>(j*(j+1)/2+i); |
865 |
|
|
} |
866 |
|
|
|
867 |
|
✗ |
return id; |
868 |
|
|
} |
869 |
|
|
|
870 |
|
✗ |
int EigenProblemPOD::thirdOrderGlobalIndex(int& i, int& j, int& k) { |
871 |
|
|
int id; |
872 |
|
|
int s[3]; |
873 |
|
✗ |
s[0] = i; |
874 |
|
✗ |
s[1] = j; |
875 |
|
✗ |
s[2] = k; |
876 |
|
|
int temp; |
877 |
|
✗ |
for (int ii=0; ii<3; ii++) { |
878 |
|
✗ |
for (int jj=0; jj<2; jj++) { |
879 |
|
✗ |
if (s[jj] < s[jj+1]) { |
880 |
|
✗ |
temp = s[jj]; |
881 |
|
✗ |
s[jj] = s[jj+1]; |
882 |
|
✗ |
s[jj+1] = temp; |
883 |
|
|
} |
884 |
|
|
} |
885 |
|
|
} |
886 |
|
✗ |
id = static_cast<int>(s[0]*(s[0]+1)*(s[0]+2)/6); |
887 |
|
✗ |
id += static_cast<int>(s[1]*(s[1]+1)/2); |
888 |
|
✗ |
id += s[2]; |
889 |
|
✗ |
return id; |
890 |
|
|
} |
891 |
|
|
|
892 |
|
✗ |
int EigenProblemPOD::fourthOrderGlobalIndex(int& i, int& j, int& k, int& h) { |
893 |
|
|
int id; |
894 |
|
|
int s[4]; |
895 |
|
✗ |
s[0] = i; |
896 |
|
✗ |
s[1] = j; |
897 |
|
✗ |
s[2] = k; |
898 |
|
✗ |
s[3] = h; |
899 |
|
|
int temp; |
900 |
|
✗ |
for (int ii=0; ii<4; ii++) { |
901 |
|
✗ |
for (int jj=0; jj<3; jj++) { |
902 |
|
✗ |
if (s[jj] < s[jj+1]) { |
903 |
|
✗ |
temp = s[jj]; |
904 |
|
✗ |
s[jj] = s[jj+1]; |
905 |
|
✗ |
s[jj+1] = temp; |
906 |
|
|
} |
907 |
|
|
} |
908 |
|
|
} |
909 |
|
✗ |
id = static_cast<int>(s[0]*(s[0]+1)*(s[0]*s[0]+5*s[0]+6)/24); |
910 |
|
✗ |
id += static_cast<int>(s[1]*(s[1]+1)*(s[1]+2)/6); |
911 |
|
✗ |
id += static_cast<int>(s[2]*(s[2]+1)/2); |
912 |
|
✗ |
id += s[3]; |
913 |
|
✗ |
return id; |
914 |
|
|
} |
915 |
|
|
|
916 |
|
|
|
917 |
|
✗ |
void EigenProblemPOD::computeTensor() { |
918 |
|
✗ |
int n = m_dimRomBasis; |
919 |
|
|
|
920 |
|
✗ |
if ( !FelisceParam::instance().hasInfarct) { |
921 |
|
|
// T_ijk = < phi_i, zeta_jk > |
922 |
|
✗ |
m_size3 = static_cast<int>(n*(n+1)*(n+2)/6); |
923 |
|
✗ |
if (m_tensorT == nullptr) |
924 |
|
✗ |
m_tensorT = new double [m_size3]; |
925 |
|
|
// Y_ijkh = < zeta_ij, zeta_kh > |
926 |
|
✗ |
m_size4 = static_cast<int>(n*(n+1)*(n*n+5*n+6)/24); |
927 |
|
✗ |
if (m_tensorY == nullptr) |
928 |
|
✗ |
m_tensorY = new double [m_size4]; |
929 |
|
✗ |
if (m_tensorE == nullptr) |
930 |
|
✗ |
m_tensorE = new double [m_size2]; |
931 |
|
|
|
932 |
|
|
|
933 |
|
|
double valueT; |
934 |
|
|
double valueY; |
935 |
|
|
double valueE; |
936 |
|
✗ |
PetscVector phiG; |
937 |
|
✗ |
phiG.duplicateFrom(m_basisPOD[0]); |
938 |
|
✗ |
PetscVector phiK; |
939 |
|
✗ |
phiK.duplicateFrom(m_basisPOD[0]); |
940 |
|
|
|
941 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
942 |
|
|
// phiK = phi_i^T * K |
943 |
|
✗ |
mult(m_Matrix[m_idK],m_basisPOD[i],phiK); |
944 |
|
✗ |
for (int j=0; j<i+1; j++) { |
945 |
|
|
// phiG = zeta_ij^T * G |
946 |
|
✗ |
mult(m_Matrix[m_idG],m_zeta[secondOrderGlobalIndex(i,j)],phiG); |
947 |
|
|
// valueE |
948 |
|
✗ |
dot(m_basisPOD[j],phiK,&valueE); |
949 |
|
✗ |
m_tensorE[secondOrderGlobalIndex(i,j)] = valueE; |
950 |
|
|
|
951 |
|
✗ |
for (int k=0; k<j+1; k++) { |
952 |
|
|
// valueT = zeta_ij^T * G * phi_k |
953 |
|
✗ |
dot(m_basisPOD[k],phiG,&valueT); |
954 |
|
✗ |
m_tensorT[thirdOrderGlobalIndex(i,j,k)] = valueT; |
955 |
|
✗ |
for (int h=0; h<k+1; h++) { |
956 |
|
|
// valueY = zeta_ij^T * G * zeta_kh |
957 |
|
✗ |
dot(m_zeta[secondOrderGlobalIndex(k,h)],phiG,&valueY); |
958 |
|
✗ |
m_tensorY[fourthOrderGlobalIndex(i,j,k,h)] = valueY; |
959 |
|
|
} |
960 |
|
|
} |
961 |
|
|
} |
962 |
|
|
} |
963 |
|
✗ |
phiG.destroy(); |
964 |
|
✗ |
} else if (FelisceParam::instance().hasInfarct) { |
965 |
|
|
// Cs_ij = < s phi_i, phi_j > |
966 |
|
✗ |
if (m_tensorCs == nullptr) |
967 |
|
✗ |
m_tensorCs = new double [m_size2]; |
968 |
|
|
// T_ijk = < phi_i, zeta_jk > |
969 |
|
✗ |
m_size3 = static_cast<int>(n*(n+1)*(n+2)/6); |
970 |
|
✗ |
if (m_tensorT == nullptr) |
971 |
|
✗ |
m_tensorT = new double [m_size3]; |
972 |
|
|
// Ts_ijk = < s phi_i, zeta_jk > |
973 |
|
✗ |
m_size3 = static_cast<int>(n*(n+1)*(n+2)/6); |
974 |
|
✗ |
if (m_tensorTs == nullptr) |
975 |
|
✗ |
m_tensorTs = new double [m_size3]; |
976 |
|
|
// Y_ijkh = < zeta_ij, zeta_kh > |
977 |
|
✗ |
m_size4 = static_cast<int>(n*(n+1)*(n*n+5*n+6)/24); |
978 |
|
✗ |
if (m_tensorY == nullptr) |
979 |
|
✗ |
m_tensorY = new double [m_size4]; |
980 |
|
✗ |
if (m_tensorE == nullptr) |
981 |
|
✗ |
m_tensorE = new double [m_size2]; |
982 |
|
|
|
983 |
|
|
double valueCs; |
984 |
|
|
double valueT; |
985 |
|
|
double valueTs; |
986 |
|
|
double valueY; |
987 |
|
|
double valueE; |
988 |
|
✗ |
PetscVector phiGs; |
989 |
|
✗ |
PetscVector zetaG; |
990 |
|
✗ |
PetscVector zetaGs; |
991 |
|
✗ |
PetscVector phiK; |
992 |
|
✗ |
phiGs.duplicateFrom(m_basisPOD[0]); |
993 |
|
✗ |
zetaG.duplicateFrom(m_basisPOD[0]); |
994 |
|
✗ |
zetaGs.duplicateFrom(m_basisPOD[0]); |
995 |
|
✗ |
phiK.duplicateFrom(m_basisPOD[0]); |
996 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
997 |
|
|
// phiGs = phi_i^T * Gs |
998 |
|
✗ |
mult(m_Matrix[m_idGs],m_basisPOD[i],phiGs); |
999 |
|
✗ |
mult(m_Matrix[m_idK],m_basisPOD[i],phiK); |
1000 |
|
✗ |
for (int j=0; j<i+1; j++) { |
1001 |
|
|
// valueCs = phi_i^T * Gs * phi_j |
1002 |
|
✗ |
dot(m_basisPOD[j],phiGs,&valueCs); |
1003 |
|
✗ |
m_tensorCs[secondOrderGlobalIndex(i,j)] = valueCs; |
1004 |
|
|
|
1005 |
|
|
//valueE |
1006 |
|
✗ |
dot(m_basisPOD[j],phiK,&valueE); |
1007 |
|
✗ |
m_tensorE[secondOrderGlobalIndex(i,j)] = valueE; |
1008 |
|
|
|
1009 |
|
|
// zetaG = zeta_ij^T * G |
1010 |
|
✗ |
mult(m_Matrix[m_idG],m_zeta[secondOrderGlobalIndex(i,j)],zetaG); |
1011 |
|
|
// zetaGs = zeta_i^T * Gs |
1012 |
|
✗ |
mult(m_Matrix[m_idGs],m_zeta[secondOrderGlobalIndex(i,j)],zetaGs); |
1013 |
|
✗ |
for (int k=0; k<j+1; k++) { |
1014 |
|
|
// valueT = zeta_ij^T * G * phi_k |
1015 |
|
✗ |
dot(m_basisPOD[k],zetaG,&valueT); |
1016 |
|
✗ |
m_tensorT[thirdOrderGlobalIndex(i,j,k)] = valueT; |
1017 |
|
|
// valueTs = zeta_ij^T * Gs * phi_k |
1018 |
|
✗ |
dot(m_basisPOD[k],zetaGs,&valueTs); |
1019 |
|
✗ |
m_tensorTs[thirdOrderGlobalIndex(i,j,k)] = valueTs; |
1020 |
|
✗ |
for (int h=0; h<k+1; h++) { |
1021 |
|
|
// valueY = zeta_ij^T * Gs * zeta_kh |
1022 |
|
✗ |
dot(m_zeta[secondOrderGlobalIndex(k,h)],zetaGs,&valueY); |
1023 |
|
✗ |
m_tensorY[fourthOrderGlobalIndex(i,j,k,h)] = valueY; |
1024 |
|
|
} |
1025 |
|
|
} |
1026 |
|
|
} |
1027 |
|
|
} |
1028 |
|
✗ |
phiGs.destroy(); |
1029 |
|
✗ |
zetaG.destroy(); |
1030 |
|
✗ |
zetaGs.destroy(); |
1031 |
|
|
} |
1032 |
|
✗ |
}; |
1033 |
|
|
|
1034 |
|
|
|
1035 |
|
|
|
1036 |
|
✗ |
void EigenProblemPOD::updateBeta() { |
1037 |
|
✗ |
if (FelisceParam::verbose() > 10 ) PetscPrintf(PETSC_COMM_WORLD,"\nEigenProblemPOD::updateBeta()\n"); |
1038 |
|
|
|
1039 |
|
✗ |
double dt = FelisceParam::instance().timeStep; |
1040 |
|
✗ |
double newBeta[m_dimRomBasis]; |
1041 |
|
|
// beta^{n+1} = beta^n + dt * gamma^n |
1042 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
1043 |
|
✗ |
newBeta[i] = m_gamma[i]; |
1044 |
|
|
} |
1045 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
1046 |
|
✗ |
m_beta[i] += dt * newBeta[i]; |
1047 |
|
|
} |
1048 |
|
|
|
1049 |
|
✗ |
double newMu[m_dimRomBasis]; |
1050 |
|
|
// mu^{n+1} = mu^n + dt * eta^n |
1051 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
1052 |
|
✗ |
newMu[i] = m_eta[i]; |
1053 |
|
|
} |
1054 |
|
✗ |
for (int i=0; i<m_dimRomBasis; i++) { |
1055 |
|
✗ |
m_mu[i] += dt * newMu[i]; |
1056 |
|
|
} |
1057 |
|
|
} |
1058 |
|
|
|
1059 |
|
|
|
1060 |
|
|
|
1061 |
|
|
// Projection functions from FE space to RO one |
1062 |
|
✗ |
void EigenProblemPOD::projectOnBasis(PetscVector& vectorDof, double* vectorBasis, bool flagMass) { |
1063 |
|
✗ |
PetscVector tmpVec; |
1064 |
|
✗ |
tmpVec.duplicateFrom(vectorDof); |
1065 |
|
|
|
1066 |
|
✗ |
if (flagMass) { |
1067 |
|
✗ |
for (int iBasis=0; iBasis<m_dimRomBasis; iBasis++) { |
1068 |
|
✗ |
mult(m_Matrix[m_idG], m_basisPOD[iBasis], tmpVec); |
1069 |
|
✗ |
dot(tmpVec,vectorDof,&vectorBasis[iBasis]); |
1070 |
|
|
} |
1071 |
|
|
} else { |
1072 |
|
✗ |
for (int iBasis=0; iBasis<m_dimRomBasis; iBasis++) { |
1073 |
|
✗ |
dot( m_basisPOD[iBasis],vectorDof,&vectorBasis[iBasis]); |
1074 |
|
|
} |
1075 |
|
|
} |
1076 |
|
|
|
1077 |
|
✗ |
tmpVec.destroy(); |
1078 |
|
|
} |
1079 |
|
|
|
1080 |
|
|
// Projection functions from RO space to FE one |
1081 |
|
✗ |
void EigenProblemPOD::projectOnDof(double* vectorBasis, PetscVector& vectorDof, int size) { |
1082 |
|
✗ |
vectorDof.set(0.); |
1083 |
|
✗ |
for (int iBasis=0; iBasis<size; iBasis++) { |
1084 |
|
✗ |
vectorDof.axpy(vectorBasis[iBasis],m_basisPOD[iBasis]); |
1085 |
|
|
} |
1086 |
|
|
} |
1087 |
|
|
|
1088 |
|
|
|
1089 |
|
|
// Auxiliary function to debug purposes. |
1090 |
|
✗ |
void EigenProblemPOD::viewPOD(double* vToPrint, int vSize, std::ofstream outStream) { |
1091 |
|
|
int rankProc; |
1092 |
|
✗ |
MPI_Comm_rank(m_petscComm,&rankProc); |
1093 |
|
✗ |
if (rankProc == 0) { |
1094 |
|
✗ |
for(int i=0; i<vSize; i++) { |
1095 |
|
✗ |
outStream << vToPrint[i] << std::endl; |
1096 |
|
|
} |
1097 |
|
✗ |
outStream.close(); |
1098 |
|
|
} |
1099 |
|
|
} |
1100 |
|
|
|
1101 |
|
|
|
1102 |
|
✗ |
void EigenProblemPOD::viewPOD(double* vToPrint, int vSize, std::string fName) { |
1103 |
|
|
int rankProc; |
1104 |
|
✗ |
MPI_Comm_rank(m_petscComm,&rankProc); |
1105 |
|
✗ |
if (rankProc == 0) { |
1106 |
|
✗ |
std::ofstream outStream(fName.c_str()); |
1107 |
|
✗ |
for(int i=0; i<vSize; i++) { |
1108 |
|
✗ |
outStream << vToPrint[i] << std::endl; |
1109 |
|
|
} |
1110 |
|
✗ |
outStream.close(); |
1111 |
|
|
} |
1112 |
|
|
} |
1113 |
|
|
|
1114 |
|
|
} |
1115 |
|
|
|