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: A. Collin |
13 |
|
|
// |
14 |
|
|
|
15 |
|
|
// System includes |
16 |
|
|
|
17 |
|
|
// External includes |
18 |
|
|
|
19 |
|
|
// Project includes |
20 |
|
|
#include "Solver/linearProblemBidomainCurv.hpp" |
21 |
|
|
#include "FiniteElement/elementMatrix.hpp" |
22 |
|
|
#include "InputOutput/io.hpp" |
23 |
|
|
|
24 |
|
|
namespace felisce |
25 |
|
|
{ |
26 |
|
✗ |
LinearProblemBidomainCurv::LinearProblemBidomainCurv(): |
27 |
|
|
LinearProblem("Problem cardiac Bidomain",2), |
28 |
|
✗ |
m_sortedSol(nullptr), |
29 |
|
✗ |
allocateSeqVec(false), |
30 |
|
✗ |
allocateLevelSet(false), |
31 |
|
✗ |
allocateSol_n_1(false), |
32 |
|
✗ |
m_avHeavU0(0.), |
33 |
|
✗ |
m_avHeavMinusU0(0.) { |
34 |
|
✗ |
if (allocateSeqVec) { |
35 |
|
✗ |
m_seqIon.destroy(); |
36 |
|
✗ |
m_seqBdfRHS.destroy(); |
37 |
|
|
} |
38 |
|
✗ |
if (allocateLevelSet) { |
39 |
|
✗ |
m_levelSet.destroy(); |
40 |
|
✗ |
m_levelSetSeq.destroy(); |
41 |
|
|
} |
42 |
|
✗ |
if (allocateSol_n_1) { |
43 |
|
✗ |
m_sol_n_1.destroy(); |
44 |
|
|
} |
45 |
|
|
} |
46 |
|
|
|
47 |
|
✗ |
LinearProblemBidomainCurv::~LinearProblemBidomainCurv() { |
48 |
|
✗ |
delete [] m_sortedSol; |
49 |
|
|
} |
50 |
|
|
|
51 |
|
✗ |
void LinearProblemBidomainCurv::initialize(std::vector<GeometricMeshRegion::Pointer>& mesh, FelisceTransient::Pointer fstransient, MPI_Comm& comm, bool doUseSNES) { |
52 |
|
✗ |
LinearProblem::initialize(mesh,comm, doUseSNES); |
53 |
|
✗ |
m_fstransient = fstransient; |
54 |
|
|
|
55 |
|
✗ |
std::vector<PhysicalVariable> listVariable(2); |
56 |
|
✗ |
std::vector<std::size_t> listNumComp(2); |
57 |
|
|
|
58 |
|
✗ |
listVariable[0] = potTransMemb; |
59 |
|
✗ |
listNumComp[0] = 1; |
60 |
|
✗ |
listVariable[1] = potExtraCell; |
61 |
|
✗ |
listNumComp[1] = 1; |
62 |
|
|
|
63 |
|
|
//define unknown of the linear system. |
64 |
|
✗ |
m_listUnknown.push_back(potTransMemb); |
65 |
|
✗ |
m_listUnknown.push_back(potExtraCell); |
66 |
|
✗ |
definePhysicalVariable(listVariable,listNumComp); |
67 |
|
|
} |
68 |
|
|
|
69 |
|
✗ |
void LinearProblemBidomainCurv::readData(IO& io) { |
70 |
|
✗ |
m_vectorFiber.resize(m_mesh[m_currentMesh]->numPoints()*3); |
71 |
|
✗ |
io.readVariable(0, 0.,&m_vectorFiber[0], m_vectorFiber.size()); |
72 |
|
✗ |
m_vectorAngle.resize(m_mesh[m_currentMesh]->numPoints()); |
73 |
|
✗ |
io.readVariable(1, 0.,&m_vectorAngle[0],m_vectorAngle.size()); |
74 |
|
✗ |
if (!FelisceParam::instance().hasHeteroCondAtria) { |
75 |
|
✗ |
m_vectorRef.resize(m_mesh[m_currentMesh]->numPoints()); //m_numDof; |
76 |
|
✗ |
io.readVariable(2, 0.,&m_vectorRef[0],m_vectorRef.size()); |
77 |
|
|
} |
78 |
|
|
} |
79 |
|
|
|
80 |
|
✗ |
void LinearProblemBidomainCurv::readDataForDA(IO& io, double iteration) { |
81 |
|
✗ |
m_data.resize(m_mesh[m_currentMesh]->numPoints()); |
82 |
|
✗ |
if (FelisceParam::instance().typeOfAppliedCurrent != "atria") |
83 |
|
✗ |
io.readVariable(2, iteration*m_fstransient->timeStep, &m_data[0],m_data.size()); |
84 |
|
|
else |
85 |
|
✗ |
io.readVariable(3, iteration*m_fstransient->timeStep, &m_data[0],m_data.size()); |
86 |
|
|
} |
87 |
|
|
|
88 |
|
✗ |
void LinearProblemBidomainCurv::getFiberDirection(felInt iel, int iUnknown, std::vector<double>& elemFiber) { |
89 |
|
|
//Warning: the fibers must be normalized. |
90 |
|
✗ |
int numSupport = static_cast<int>(supportDofUnknown(iUnknown).iEle()[iel+1] - supportDofUnknown(iUnknown).iEle()[iel]); |
91 |
|
✗ |
elemFiber.resize( numSupport*3,0.); |
92 |
|
✗ |
for (int i = 0; i < numSupport; i++) { |
93 |
|
✗ |
elemFiber[i*3] = m_vectorFiber[3*(supportDofUnknown(iUnknown).iSupportDof()[supportDofUnknown(iUnknown).iEle()[iel]+i])]; |
94 |
|
✗ |
elemFiber[i*3+1] = m_vectorFiber[3*(supportDofUnknown(iUnknown).iSupportDof()[supportDofUnknown(iUnknown).iEle()[iel]+i])+1]; |
95 |
|
✗ |
elemFiber[i*3+2] = m_vectorFiber[3*(supportDofUnknown(iUnknown).iSupportDof()[supportDofUnknown(iUnknown).iEle()[iel]+i])+2]; |
96 |
|
|
} |
97 |
|
|
} |
98 |
|
|
|
99 |
|
✗ |
void LinearProblemBidomainCurv::getAngle(felInt iel, int iUnknown, std::vector<double>& elemAngle) { |
100 |
|
✗ |
int numSupport = static_cast<int>(supportDofUnknown(iUnknown).iEle()[iel+1] - supportDofUnknown(iUnknown).iEle()[iel]); |
101 |
|
✗ |
elemAngle.resize(numSupport,0.); |
102 |
|
✗ |
for (int i = 0; i < numSupport; i++) { |
103 |
|
✗ |
elemAngle[i] = m_vectorAngle[(supportDofUnknown(iUnknown).iSupportDof()[supportDofUnknown(iUnknown).iEle()[iel]+i])]; |
104 |
|
|
} |
105 |
|
|
} |
106 |
|
|
|
107 |
|
✗ |
void LinearProblemBidomainCurv::getData(felInt iel, int iUnknown, std::vector<double>& elemData) { |
108 |
|
|
//double& vGate = FelisceParam::instance().vGate; |
109 |
|
✗ |
double vGate = -67.; |
110 |
|
✗ |
int numSupport = static_cast<int>(supportDofUnknown(iUnknown).iEle()[iel+1] - supportDofUnknown(iUnknown).iEle()[iel]); |
111 |
|
✗ |
elemData.clear(); |
112 |
|
✗ |
elemData.resize(numSupport,0.); |
113 |
|
✗ |
for (int i = 0; i < numSupport; i++) { |
114 |
|
✗ |
if (m_data[(supportDofUnknown(iUnknown).iSupportDof()[supportDofUnknown(iUnknown).iEle()[iel]+i])] > vGate) { |
115 |
|
✗ |
elemData[i] = 1.0; |
116 |
|
|
} else |
117 |
|
✗ |
elemData[i] = -1.0; |
118 |
|
|
} |
119 |
|
|
} |
120 |
|
|
|
121 |
|
✗ |
void LinearProblemBidomainCurv::initPerElementTypeBD(ElementType eltType, FlagMatrixRHS flagMatrixRHS) { |
122 |
|
|
IGNORE_UNUSED_FLAG_MATRIX_RHS; |
123 |
|
|
IGNORE_UNUSED_ELT_TYPE; |
124 |
|
|
//Init pointer on Finite Element, Variable or idVariable |
125 |
|
✗ |
m_ipotTransMemb = m_listVariable.getVariableIdList(potTransMemb); |
126 |
|
✗ |
m_ipotExtraCell = m_listVariable.getVariableIdList(potExtraCell); |
127 |
|
✗ |
m_fePotTransMemb = m_listCurvilinearFiniteElement[m_ipotTransMemb]; |
128 |
|
✗ |
m_fePotExtraCell = m_fePotTransMemb; |
129 |
|
✗ |
if (FelisceParam::instance().dataAssimilation) { |
130 |
|
✗ |
m_elemFieldVm.initialize(DOF_FIELD,*m_fePotTransMemb,1); |
131 |
|
✗ |
m_elemFieldLSVm.initialize(DOF_FIELD,*m_fePotTransMemb,1); |
132 |
|
✗ |
m_avHeavU0 = averageLS(m_fePotTransMemb,m_ipotTransMemb,0); |
133 |
|
✗ |
m_avHeavMinusU0 = averageLS(m_fePotTransMemb,m_ipotTransMemb,1); |
134 |
|
|
} |
135 |
|
|
} |
136 |
|
|
|
137 |
|
✗ |
double LinearProblemBidomainCurv::conductivityHom(int currentLabel) { |
138 |
|
✗ |
double conductivity = 1.0; |
139 |
|
|
|
140 |
|
✗ |
if (FelisceParam::instance().hasHeteroCondAtria) { |
141 |
|
✗ |
if (std::fabs(currentLabel - 31.0) < 1.0e-12) { //SN |
142 |
|
✗ |
conductivity = 0.3; |
143 |
|
✗ |
} else if (std::fabs(currentLabel - 32.0) < 1.0e-12) { //CT |
144 |
|
✗ |
conductivity = 3.0; |
145 |
|
✗ |
} else if (std::fabs(currentLabel - 33.0) < 1.0e-12) { //BB |
146 |
|
✗ |
conductivity = 4.75; |
147 |
|
✗ |
} else if (std::fabs(currentLabel - 34.0) < 1.0e-12) { //RAI |
148 |
|
✗ |
conductivity = 0.375; |
149 |
|
✗ |
} else if (std::fabs(currentLabel - 35.0) < 1.0e-12) { //PM |
150 |
|
✗ |
conductivity = 1.5; |
151 |
|
|
} else { //regular atria and F0 |
152 |
|
✗ |
conductivity = 1.0; |
153 |
|
|
} |
154 |
|
|
} |
155 |
|
|
|
156 |
|
✗ |
return conductivity; |
157 |
|
|
} |
158 |
|
|
|
159 |
|
✗ |
double LinearProblemBidomainCurv::conductivityHet(int currentLabel) { |
160 |
|
✗ |
double conductivity = 1.0; |
161 |
|
|
|
162 |
|
✗ |
if (FelisceParam::instance().hasHeteroCondAtria) { |
163 |
|
✗ |
if (std::fabs(currentLabel - 31.0) < 1.0e-12) { //SN |
164 |
|
✗ |
conductivity = 0.36; |
165 |
|
✗ |
} else if (std::fabs(currentLabel - 32.0) < 1.0e-12) { //CT |
166 |
|
✗ |
conductivity = 4.5; |
167 |
|
✗ |
} else if (std::fabs(currentLabel - 33.0) < 1.0e-12) { //BB |
168 |
|
✗ |
conductivity = 7.75; |
169 |
|
✗ |
} else if (std::fabs(currentLabel - 34.0) < 1.0e-12) { //RAI |
170 |
|
✗ |
conductivity = 0.45; |
171 |
|
✗ |
} else if (std::fabs(currentLabel - 35.0) < 1.0e-12) { //PM |
172 |
|
✗ |
conductivity = 1.8; |
173 |
|
✗ |
} else if (std::fabs(currentLabel - 36.0) < 1.0e-12) { //F0 |
174 |
|
✗ |
conductivity = 0.9; |
175 |
|
|
} else { //regular atria |
176 |
|
✗ |
conductivity = 1.0; |
177 |
|
|
} |
178 |
|
|
} |
179 |
|
|
|
180 |
|
✗ |
return conductivity; |
181 |
|
|
} |
182 |
|
|
|
183 |
|
✗ |
void LinearProblemBidomainCurv::initLevelSet() { |
184 |
|
✗ |
m_levelSet.duplicateFrom(solution()); |
185 |
|
|
|
186 |
|
✗ |
m_levelSet.zeroEntries(); |
187 |
|
✗ |
m_levelSet.assembly(); |
188 |
|
|
|
189 |
|
✗ |
m_levelSetSeq.createSeq(PETSC_COMM_SELF,m_numDof); |
190 |
|
✗ |
m_levelSetSeq.set(0.0); |
191 |
|
|
|
192 |
|
✗ |
allocateLevelSet = true; |
193 |
|
|
} |
194 |
|
|
|
195 |
|
✗ |
void LinearProblemBidomainCurv::levelSet(double min, double max) { |
196 |
|
|
//double& vGate = FelisceParam::instance().vGate; |
197 |
|
✗ |
double vGate = -67.; |
198 |
|
|
|
199 |
|
|
felInt pos; |
200 |
|
|
double value_sol; |
201 |
|
|
double value_levelSet; |
202 |
|
|
|
203 |
|
✗ |
for (felInt i = 0; i < numDofLocalPerUnknown(potTransMemb) ; i++) { |
204 |
|
✗ |
ISLocalToGlobalMappingApply(mappingDofLocalToDofGlobal(potTransMemb),1,&i,&pos); |
205 |
|
✗ |
solution().getValues(1,&pos,&value_sol); |
206 |
|
✗ |
if( value_sol > vGate) { |
207 |
|
✗ |
value_levelSet = max; |
208 |
|
✗ |
} else if( value_sol < vGate) { |
209 |
|
✗ |
value_levelSet = min; |
210 |
|
|
} else |
211 |
|
✗ |
value_levelSet = 0.0; |
212 |
|
|
|
213 |
|
✗ |
m_levelSet.setValue(pos,value_levelSet, INSERT_VALUES); |
214 |
|
|
} |
215 |
|
✗ |
m_levelSet.assembly(); |
216 |
|
|
|
217 |
|
✗ |
gatherVector(m_levelSet, m_levelSetSeq); |
218 |
|
|
} |
219 |
|
|
|
220 |
|
✗ |
void LinearProblemBidomainCurv::initSolutionTimeBef() { |
221 |
|
✗ |
m_sol_n_1.duplicateFrom(solution()); |
222 |
|
✗ |
m_sol_n_1.copyFrom(solution()); |
223 |
|
|
|
224 |
|
✗ |
allocateSol_n_1 = true; |
225 |
|
|
} |
226 |
|
|
|
227 |
|
✗ |
void LinearProblemBidomainCurv::solutionTimeBef() { |
228 |
|
✗ |
m_sol_n_1.copyFrom(solution()); |
229 |
|
|
} |
230 |
|
|
|
231 |
|
✗ |
double LinearProblemBidomainCurv::averageLS(CurvilinearFiniteElement* fePtr, int idfe, int minus) { |
232 |
|
✗ |
FEL_ASSERT(fePtr) |
233 |
|
✗ |
CurvilinearFiniteElement& fe = *fePtr; |
234 |
|
|
|
235 |
|
✗ |
Variable variablePotTransMemb = m_listVariable[idfe]; |
236 |
|
|
|
237 |
|
✗ |
int numComp = variablePotTransMemb.numComponent(); |
238 |
|
|
// dofValue contains : values in dof by element. |
239 |
|
✗ |
double dofValueSeqSol[fe.numDof()*numComp]; |
240 |
|
|
|
241 |
|
|
// idGlobalDof contains: global number od dof by element. |
242 |
|
✗ |
felInt idGlobalDof[fe.numDof()*numComp]; |
243 |
|
|
|
244 |
|
|
//geometric element type in the mesh. |
245 |
|
|
ElementType eltType; |
246 |
|
|
//number of points per geometric element. |
247 |
|
✗ |
int numPointPerElt = 0; |
248 |
|
|
//number of element for one label and one eltType. |
249 |
|
✗ |
felInt numEltPerLabel = 0; |
250 |
|
|
//Points of the current element. |
251 |
|
✗ |
std::vector<Point*> elemPoint; |
252 |
|
|
//Id of points of the element. |
253 |
|
✗ |
std::vector<felInt> elemIdPoint; |
254 |
|
|
|
255 |
|
|
// use to identify global id of dof. |
256 |
|
|
felInt idDof; |
257 |
|
✗ |
felInt cpt = 0; |
258 |
|
✗ |
double integralValue = 0.; |
259 |
|
✗ |
double heavisideMeasure = 0.; |
260 |
|
✗ |
double elementValueHeav = 0.; |
261 |
|
✗ |
double elementValueHeavU0 = 0.; |
262 |
|
|
//Id link element to its supportDof. |
263 |
|
|
felInt ielSupportDof; |
264 |
|
|
|
265 |
|
✗ |
double vGate = -67.; |
266 |
|
|
|
267 |
|
|
// initialize global numbering from mesh. |
268 |
|
|
felInt numElement[ GeometricMeshRegion::m_numTypesOfElement ]; |
269 |
|
✗ |
for (int ityp=0; ityp<GeometricMeshRegion::m_numTypesOfElement; ityp++ ) { |
270 |
|
✗ |
eltType = (ElementType)ityp; |
271 |
|
✗ |
numElement[eltType] = 0; |
272 |
|
|
} |
273 |
|
|
|
274 |
|
|
// Loop on elementType. |
275 |
|
✗ |
const std::vector<ElementType>& bagElementTypeDomain = m_meshLocal[m_currentMesh]->bagElementTypeDomain(); |
276 |
|
✗ |
for (std::size_t i = 0; i < bagElementTypeDomain.size(); ++i) { |
277 |
|
✗ |
eltType = bagElementTypeDomain[i]; |
278 |
|
✗ |
numPointPerElt = GeometricMeshRegion::m_numPointsPerElt[eltType]; |
279 |
|
✗ |
elemPoint.resize(numPointPerElt, nullptr); |
280 |
|
✗ |
elemIdPoint.resize(numPointPerElt,0); |
281 |
|
|
// Loop on region define in the mesh. |
282 |
|
|
|
283 |
|
✗ |
for(auto itRef = m_meshLocal[m_currentMesh]->intRefToBegEndMaps[eltType].begin(); |
284 |
|
✗ |
itRef != m_meshLocal[m_currentMesh]->intRefToBegEndMaps[eltType].end(); itRef++) { |
285 |
|
✗ |
numEltPerLabel = itRef->second.second; |
286 |
|
|
|
287 |
|
|
// Loop on elements. |
288 |
|
✗ |
for ( felInt iel = 0; iel < numEltPerLabel; iel++) { |
289 |
|
|
// get points of element. |
290 |
|
✗ |
setElemPoint(eltType, numElement[eltType], elemPoint, elemIdPoint, &ielSupportDof); |
291 |
|
|
// compute measure of the current element. |
292 |
|
✗ |
fe.updateMeas(0, elemPoint); |
293 |
|
|
|
294 |
|
✗ |
cpt = 0; |
295 |
|
|
// identify global id of Dof for the current element. |
296 |
|
✗ |
for(int iComp = 0; iComp<numComp; iComp++) { |
297 |
|
✗ |
for(int iDof=0; iDof<fe.numDof(); iDof++) { |
298 |
|
✗ |
dof().loc2glob(ielSupportDof,iDof,m_ipotTransMemb,iComp,idDof); |
299 |
|
✗ |
idGlobalDof[cpt] = idDof; |
300 |
|
✗ |
cpt++; |
301 |
|
|
} |
302 |
|
|
} |
303 |
|
|
|
304 |
|
✗ |
std::vector<double> dofValueData; |
305 |
|
✗ |
getData(ielSupportDof, m_ipotTransMemb, dofValueData); |
306 |
|
|
|
307 |
|
|
//mapping to obtain new global numbering. (I/O array: idGlobalDof). |
308 |
|
✗ |
AOApplicationToPetsc(m_ao,fe.numDof()*numComp,idGlobalDof); |
309 |
|
|
//gets values of associate dofs. |
310 |
|
|
|
311 |
|
✗ |
this->sequentialSolution().getValues(fe.numDof()*numComp,idGlobalDof,dofValueSeqSol); |
312 |
|
|
|
313 |
|
|
//Interpolation |
314 |
|
✗ |
for(int ig=0; ig<fe.numQuadraturePoint(); ig++) { |
315 |
|
✗ |
elementValueHeavU0 = 0.; |
316 |
|
✗ |
elementValueHeav = 0.; |
317 |
|
✗ |
cpt = 0; |
318 |
|
✗ |
for(int icomp=0; icomp<numComp; icomp++) { |
319 |
|
✗ |
for(int j=0; j<fe.numDof(); j++) { |
320 |
|
✗ |
if (minus == 1) { |
321 |
|
✗ |
if (dofValueSeqSol[cpt] < vGate) { |
322 |
|
✗ |
elementValueHeavU0 += fe.phi[ig](j) * 2. * dofValueData[cpt]; |
323 |
|
✗ |
elementValueHeav += fe.phi[ig](j) * 2.; |
324 |
|
|
} |
325 |
|
|
} else { |
326 |
|
✗ |
if (dofValueSeqSol[cpt] > vGate) { |
327 |
|
✗ |
elementValueHeavU0 += fe.phi[ig](j) * 2. * dofValueData[cpt]; |
328 |
|
✗ |
elementValueHeav += fe.phi[ig](j) * 2.; |
329 |
|
|
} |
330 |
|
|
} |
331 |
|
|
|
332 |
|
✗ |
cpt++; |
333 |
|
|
} |
334 |
|
✗ |
integralValue += elementValueHeavU0 * fe.weightMeas(ig); |
335 |
|
✗ |
heavisideMeasure += elementValueHeav * fe.weightMeas(ig); |
336 |
|
|
} |
337 |
|
|
} |
338 |
|
✗ |
numElement[eltType]++; |
339 |
|
|
} |
340 |
|
|
} |
341 |
|
|
} |
342 |
|
|
|
343 |
|
✗ |
double valueReceive = 0.; |
344 |
|
✗ |
MPI_Allreduce(&integralValue,&valueReceive, 1,MPI_DOUBLE,MPI_SUM,MpiInfo::petscComm()); |
345 |
|
✗ |
integralValue = valueReceive; |
346 |
|
✗ |
valueReceive = 0.; |
347 |
|
✗ |
MPI_Allreduce(&heavisideMeasure,&valueReceive, 1,MPI_DOUBLE,MPI_SUM,MpiInfo::petscComm()); |
348 |
|
✗ |
heavisideMeasure = valueReceive; |
349 |
|
|
|
350 |
|
✗ |
return integralValue/(heavisideMeasure+0.000000001); |
351 |
|
|
} |
352 |
|
|
|
353 |
|
|
|
354 |
|
✗ |
void LinearProblemBidomainCurv::computeElementArrayBD(const std::vector<Point*>& elemPoint, const std::vector<felInt>& elemIdPoint, felInt& iel, FlagMatrixRHS flagMatrixRHS) { |
355 |
|
|
IGNORE_UNUSED_ELEM_ID_POINT; |
356 |
|
|
|
357 |
|
✗ |
m_fePotTransMemb->updateMeasNormalContra(0, elemPoint); |
358 |
|
|
|
359 |
|
✗ |
if (flagMatrixRHS == FlagMatrixRHS::only_rhs) { |
360 |
|
✗ |
if (FelisceParam::instance().dataAssimilation) { |
361 |
|
✗ |
std::vector<double> elemData; |
362 |
|
✗ |
getData(iel,m_ipotTransMemb,elemData); |
363 |
|
|
|
364 |
|
✗ |
m_elemFieldVm.setValue(this->sequentialSolution(), *m_fePotTransMemb, iel, m_ipotTransMemb, m_ao, dof()); |
365 |
|
|
|
366 |
|
✗ |
m_elemFieldLSVm.setValue(m_levelSetSeq, *m_fePotTransMemb, iel, m_ipotTransMemb, m_ao, dof()); |
367 |
|
|
|
368 |
|
✗ |
m_elementVectorBD[0]->grad_u_dot_grad_LSu_v(1.0, m_elemFieldVm, m_elemFieldLSVm, elemData, m_avHeavU0, m_avHeavMinusU0, FelisceParam::instance().insidePar, FelisceParam::instance().outsidePar, *m_fePotTransMemb,0); |
369 |
|
|
} |
370 |
|
|
} else { |
371 |
|
|
// Assembling matrix(0) |
372 |
|
✗ |
const double coefTime = 1./m_fstransient->timeStep; |
373 |
|
✗ |
const double coefAm = FelisceParam::instance().Am; |
374 |
|
✗ |
const double coefCm = FelisceParam::instance().Cm; |
375 |
|
✗ |
const double coef = m_bdf->coeffDeriv0()*coefTime*coefAm*coefCm; |
376 |
|
|
|
377 |
|
|
//1st equation : Am*Cm/dt V_m |
378 |
|
|
|
379 |
|
✗ |
m_elementMatBD[0]->phi_i_phi_j(coef,*m_fePotTransMemb,0,0,1); |
380 |
|
|
|
381 |
|
✗ |
if (FelisceParam::instance().testCase == 1) { |
382 |
|
✗ |
std::vector<double> elemFiber; |
383 |
|
✗ |
getFiberDirection(iel,m_ipotTransMemb,elemFiber); |
384 |
|
|
|
385 |
|
✗ |
std::vector<double> elemAngle; |
386 |
|
✗ |
getAngle(iel,m_ipotTransMemb,elemAngle); |
387 |
|
|
|
388 |
|
|
//1st equation : \sigma_i^t \grad V_m |
389 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor*conductivityHom(m_currentLabel),*m_fePotTransMemb,0,0,1); |
390 |
|
✗ |
m_elementMatBD[0]->tau_orthotau_grad_phi_i_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor)*conductivityHet(m_currentLabel),elemFiber,elemAngle,*m_fePotTransMemb,0,0,1); |
391 |
|
|
//m_elementMatBD[0]->tau_grad_phi_i_tau_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor),elemFiber,*m_fePotTransMemb,0,0,1); |
392 |
|
|
|
393 |
|
|
//1st equation : \sigma_i^t \grad u_e |
394 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor*conductivityHom(m_currentLabel),*m_fePotExtraCell,0,1,1); |
395 |
|
✗ |
m_elementMatBD[0]->tau_orthotau_grad_phi_i_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor)*conductivityHet(m_currentLabel),elemFiber,elemAngle,*m_fePotExtraCell,0,1,1); |
396 |
|
|
//m_elementMatBD[0]->tau_grad_phi_i_tau_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor),elemFiber,*m_fePotExtraCell,0,1,1); |
397 |
|
|
|
398 |
|
|
//2nd equation : \sigma_i^t \grad V_m |
399 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor*conductivityHom(m_currentLabel),*m_fePotTransMemb,1,0,1); |
400 |
|
✗ |
m_elementMatBD[0]->tau_orthotau_grad_phi_i_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor)*conductivityHet(m_currentLabel),elemFiber,elemAngle,*m_fePotTransMemb,1,0,1); |
401 |
|
|
//m_elementMatBD[0]->tau_grad_phi_i_tau_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor),elemFiber,*m_fePotTransMemb,1,0,1); |
402 |
|
|
|
403 |
|
|
//2nd equation : \eps u_e |
404 |
|
✗ |
double eps = 1.0e-6; |
405 |
|
✗ |
m_elementMatBD[0]->phi_i_phi_j(eps,*m_fePotExtraCell,1,1,1); |
406 |
|
|
|
407 |
|
|
//2nd equation : (\sigma_i^t+\sigma_e^t) \grad u_e |
408 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j((FelisceParam::instance().intraTransvTensor+FelisceParam::instance().extraTransvTensor)*conductivityHom(m_currentLabel),*m_fePotExtraCell,1,1,1); |
409 |
|
✗ |
m_elementMatBD[0]->tau_orthotau_grad_phi_i_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor + FelisceParam::instance().extraFiberTensor - FelisceParam::instance().extraTransvTensor)*conductivityHet(m_currentLabel),elemFiber,elemAngle,*m_fePotExtraCell,1,1,1); |
410 |
|
|
//m_elementMatBD[0]->tau_grad_phi_i_tau_grad_phi_j((FelisceParam::instance().intraFiberTensor - FelisceParam::instance().intraTransvTensor + FelisceParam::instance().extraFiberTensor - FelisceParam::instance().extraTransvTensor),elemFiber,*m_fePotExtraCell,1,1,1); |
411 |
|
✗ |
} else { |
412 |
|
|
//1st equation : \sigma_i^t \grad V_m |
413 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor,*m_fePotTransMemb,0,0,1); |
414 |
|
|
|
415 |
|
|
//1st equation : \sigma_i^t \grad u_e |
416 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor,*m_fePotExtraCell,0,1,1); |
417 |
|
|
|
418 |
|
|
//2nd equation : \sigma_i^t \grad V_m |
419 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j(FelisceParam::instance().intraTransvTensor,*m_fePotTransMemb,1,0,1); |
420 |
|
|
|
421 |
|
|
//2nd equation : \eps u_e |
422 |
|
✗ |
double eps = 1.0e-6; |
423 |
|
✗ |
m_elementMatBD[0]->phi_i_phi_j(eps,*m_fePotExtraCell,1,1,1); |
424 |
|
|
|
425 |
|
|
//2nd equation : (\sigma_i^t+\sigma_e^t) \grad u_e |
426 |
|
✗ |
m_elementMatBD[0]->grad_phi_i_grad_phi_j((FelisceParam::instance().intraTransvTensor+FelisceParam::instance().extraTransvTensor),*m_fePotExtraCell,1,1,1); |
427 |
|
|
|
428 |
|
|
} |
429 |
|
|
|
430 |
|
|
//Assembling matrix(1) = mass (only first quadrant if hasAppliedExteriorCurrent = false) |
431 |
|
✗ |
m_elementMatBD[1]->phi_i_phi_j(1.,*m_fePotTransMemb,0,0,1); |
432 |
|
✗ |
if (FelisceParam::instance().hasAppliedExteriorCurrent) |
433 |
|
✗ |
m_elementMatBD[1]->phi_i_phi_j(1.,*m_fePotTransMemb,1,1,1); |
434 |
|
|
|
435 |
|
✗ |
if (FelisceParam::instance().dataAssimilation) { |
436 |
|
✗ |
std::vector<double> elemData; |
437 |
|
✗ |
getData(iel,m_ipotTransMemb,elemData); |
438 |
|
✗ |
m_elemFieldVm.setValue(this->sequentialSolution(), *m_fePotTransMemb, iel, m_ipotTransMemb, m_ao, dof()); |
439 |
|
|
|
440 |
|
✗ |
m_elemFieldLSVm.setValue(m_levelSetSeq, *m_fePotTransMemb, iel, m_ipotTransMemb, m_ao, dof()); |
441 |
|
✗ |
m_elementVectorBD[0]->grad_u_dot_grad_LSu_v(1.0, m_elemFieldVm, m_elemFieldLSVm, elemData, m_avHeavU0, m_avHeavMinusU0, FelisceParam::instance().insidePar, FelisceParam::instance().outsidePar, *m_fePotTransMemb,0); |
442 |
|
|
} |
443 |
|
|
} |
444 |
|
|
|
445 |
|
|
|
446 |
|
|
} |
447 |
|
|
|
448 |
|
✗ |
void LinearProblemBidomainCurv::addMatrixRHS() { |
449 |
|
|
// _RHS = m_mass * _RHS. |
450 |
|
✗ |
PetscVector tmpB; |
451 |
|
✗ |
tmpB.duplicateFrom(vector()); |
452 |
|
✗ |
tmpB.copyFrom(vector()); |
453 |
|
✗ |
mult(matrix(1),tmpB,vector()); |
454 |
|
✗ |
tmpB.destroy(); |
455 |
|
|
} |
456 |
|
|
|
457 |
|
✗ |
void LinearProblemBidomainCurv::writeEnsightScalar(double* solValue, int idIter, std::string varName) { |
458 |
|
✗ |
std::string iteration; |
459 |
|
✗ |
if (idIter < 10) |
460 |
|
✗ |
iteration = "0000" + std::to_string(idIter); |
461 |
|
✗ |
else if (idIter < 100) |
462 |
|
✗ |
iteration = "000" + std::to_string(idIter); |
463 |
|
✗ |
else if (idIter < 1000) |
464 |
|
✗ |
iteration = "00" + std::to_string(idIter); |
465 |
|
✗ |
else if (idIter < 10000) |
466 |
|
✗ |
iteration = "0" + std::to_string(idIter); |
467 |
|
✗ |
else if (idIter < 100000) |
468 |
|
✗ |
iteration = std::to_string(idIter); |
469 |
|
|
|
470 |
|
✗ |
std::string fileName = FelisceParam::instance().resultDir + "/" + varName + "." + iteration + ".scl"; |
471 |
|
|
FILE * solFile; |
472 |
|
✗ |
solFile = fopen(fileName.c_str(),"w"); |
473 |
|
✗ |
fprintf(solFile, "Scalar per node\n"); |
474 |
|
✗ |
int count = 0; |
475 |
|
✗ |
for(int j=0; j < static_cast<int>(m_numDof/2); j++) { |
476 |
|
✗ |
fprintf(solFile,"%12.5e", solValue[j]); |
477 |
|
✗ |
count++; |
478 |
|
✗ |
if(count == 6) { |
479 |
|
✗ |
fprintf(solFile,"\n"); |
480 |
|
✗ |
count = 0; |
481 |
|
|
} |
482 |
|
|
} |
483 |
|
✗ |
fprintf(solFile,"\n"); |
484 |
|
✗ |
fclose(solFile); |
485 |
|
|
} |
486 |
|
|
|
487 |
|
✗ |
void LinearProblemBidomainCurv::writeEnsightCase(int numIt, std::string varName) { |
488 |
|
|
|
489 |
|
|
FILE * pFile; |
490 |
|
✗ |
std::string caseName = FelisceParam::instance().resultDir + "/" + varName + ".case"; |
491 |
|
✗ |
pFile = fopen(caseName.c_str(),"w"); |
492 |
|
✗ |
fprintf(pFile,"FORMAT\n"); |
493 |
|
✗ |
fprintf(pFile,"type: ensight\n"); |
494 |
|
✗ |
fprintf(pFile, "GEOMETRY\n"); |
495 |
|
✗ |
std::string nameGeo = "model: 1 " + FelisceParam::instance().outputMesh[0] + "\n"; |
496 |
|
✗ |
fprintf(pFile,"%s", nameGeo.c_str()); |
497 |
|
✗ |
fprintf(pFile, "VARIABLE\n"); |
498 |
|
✗ |
std::string secName = "scalar per node: 1 " + varName + " " + varName +".*****.scl\n"; |
499 |
|
✗ |
fprintf(pFile, "%s", secName.c_str()); |
500 |
|
✗ |
fprintf(pFile, "TIME\n"); |
501 |
|
✗ |
fprintf(pFile, "time std::set: %d \n", 1); |
502 |
|
✗ |
fprintf(pFile, "number of steps: %d \n", numIt); |
503 |
|
✗ |
fprintf(pFile, "filename start number: %d\n", 0); |
504 |
|
✗ |
fprintf(pFile, "filename increment: %d\n", 1); |
505 |
|
✗ |
fprintf(pFile, "time values:\n"); |
506 |
|
|
|
507 |
|
✗ |
const double dt = FelisceParam::instance().timeStep * FelisceParam::instance().frequencyWriteSolution; |
508 |
|
✗ |
int count=0; |
509 |
|
✗ |
for(int j=0; j < numIt; j++) { |
510 |
|
✗ |
double valT = j*dt; |
511 |
|
✗ |
fprintf(pFile,"%lf ",valT); |
512 |
|
✗ |
count++; |
513 |
|
✗ |
if(count == 6) { |
514 |
|
✗ |
fprintf(pFile,"\n"); |
515 |
|
✗ |
count=0; |
516 |
|
|
} |
517 |
|
|
} |
518 |
|
✗ |
fclose(pFile); |
519 |
|
|
} |
520 |
|
|
} |
521 |
|
|
|