GCC Code Coverage Report


Directory: ./
File: Model/NSSimplifiedFSIModel.cpp
Date: 2024-04-14 07:32:34
Exec Total Coverage
Lines: 141 156 90.4%
Branches: 103 194 53.1%

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:
13 //
14
15 // System includes
16
17 // External includes
18
19 // Project includes
20 #include "Model/NSSimplifiedFSIModel.hpp"
21 #include "Solver/linearProblemNSSimplifiedFSI.hpp"
22
23 namespace felisce {
24
25
4/8
✓ Branch 2 taken 12 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 12 times.
✗ Branch 6 not taken.
✓ Branch 10 taken 12 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 12 times.
✗ Branch 14 not taken.
12 NSSimplifiedFSIModel::NSSimplifiedFSIModel():Model() {
26 12 m_lumpedModelBC = nullptr;
27
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 m_name = "NS";
28 12 m_exportBdQuantities=true;
29 12 m_idNSSimpl = 0; // in this class is obviously zero, but what if you have derived model?
30 12 m_SimplifiedFSI=nullptr;
31 12 }
32
33
34 24 NSSimplifiedFSIModel::~NSSimplifiedFSIModel() {
35
1/2
✓ Branch 0 taken 12 times.
✗ Branch 1 not taken.
24 if( m_SimplifiedFSI )
36
2/2
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 8 times.
24 if ( m_SimplifiedFSI->autoregulated() ) {
37 8 m_SimplifiedFSI->exportControl(MpiInfo::rankProc());
38 }
39
2/2
✓ Branch 0 taken 4 times.
✓ Branch 1 taken 8 times.
24 if(m_lumpedModelBC)
40
1/2
✓ Branch 0 taken 4 times.
✗ Branch 1 not taken.
8 delete m_lumpedModelBC;
41 }
42
43
44 12 void NSSimplifiedFSIModel::initializeDerivedModel() {
45
2/2
✓ Branch 2 taken 4 times.
✓ Branch 3 taken 8 times.
12 if( FelisceParam::instance().lumpedModelBCLabel.size() ) {
46
2/4
✓ Branch 2 taken 4 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 4 times.
✗ Branch 6 not taken.
4 m_lumpedModelBC = new LumpedModelBC( m_fstransient );
47 4 static_cast<LinearProblemNSSimplifiedFSI*>(m_linearProblem[m_idNSSimpl])->initializeLumpedModelBC(m_lumpedModelBC);
48 }
49
50
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 12 times.
12 if (FelisceParam::instance().orderBdfNS > 2)
51 FEL_ERROR("BDF not yet implemented for order greater than 2 with Navier-Stokes.");
52
53 12 m_bdfNS.defineOrder(FelisceParam::instance().orderBdfNS);
54 12 m_linearProblem[m_idNSSimpl]->initializeTimeScheme(&m_bdfNS);
55 12 }
56
57
58 16 void NSSimplifiedFSIModel::preAssemblingMatrixRHS(std::size_t iProblem) {
59
2/2
✓ Branch 0 taken 12 times.
✓ Branch 1 taken 4 times.
16 if ( iProblem == m_idNSSimpl) {
60 12 m_solExtrapolate.duplicateFrom(m_linearProblem[iProblem]->vector());
61 12 m_solExtrapolate.zeroEntries();
62 }
63 16 m_linearProblem[iProblem]->solution().zeroEntries();
64 16 m_linearProblem[iProblem]->assembleMatrixRHS(MpiInfo::rankProc());
65 16 }
66
67 64 void NSSimplifiedFSIModel::postAssemblingMatrixRHS(std::size_t iProblem) {
68 // _Matrix = _A + _Matrix (add static matrix to the dynamic matrix to build
69 // complete matrix of the system
70 64 m_linearProblem[iProblem]->addMatrixRHS();
71 64 }
72
73 12 void NSSimplifiedFSIModel::writeExtraVectors() {
74 12 m_SimplifiedFSI = static_cast<LinearProblemNSSimplifiedFSI*>(m_linearProblem[m_idNSSimpl]);
75
1/2
✓ Branch 0 taken 12 times.
✗ Branch 1 not taken.
12 if ( m_exportBdQuantities ) {
76
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 PetscVector zero;
77
1/2
✓ Branch 3 taken 12 times.
✗ Branch 4 not taken.
12 zero.duplicateFrom(m_linearProblem[m_idNSSimpl]->solution());
78
1/2
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
12 zero.zeroEntries();
79
3/6
✓ Branch 3 taken 12 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 12 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 12 times.
✗ Branch 12 not taken.
12 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("displacement"));
80
81
3/6
✓ Branch 3 taken 12 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 12 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 12 times.
✗ Branch 12 not taken.
12 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("normal"));
82
2/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 12 times.
12 if ( FelisceParam::instance().exportDofPartition )
83 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("partition"));
84
85
3/4
✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 4 times.
12 if (FelisceParam::instance().exportPrincipalDirectionsAndCurvature ) {
86
3/6
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 8 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
8 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("maxEigenVector"));
87
3/6
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 8 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
8 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("minEigenVector"));
88
3/6
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 8 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
8 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("meanCurvature"));
89
3/6
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 8 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
8 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec( zero , MpiInfo::rankProc(), m_ios, m_fstransient->time, m_fstransient->iteration, std::string("koiterCoeff"));
90 }
91 12 }
92 12 }
93
94
95 //Everything in this function should depend only on the previous time step and NOT on current data!!
96 64 void NSSimplifiedFSIModel::prepareForward(){
97
2/2
✓ Branch 1 taken 12 times.
✓ Branch 2 taken 52 times.
64 if( m_fstransient->iteration == 0 ) {
98 12 writeExtraVectors();
99 }
100
101 // Write solution for postprocessing (if required)
102 64 writeSolution();
103
104 // Advance time step.
105 64 updateTime();
106
107 // Print time information
108 64 printNewTimeIterationBanner();
109
110
2/2
✓ Branch 1 taken 12 times.
✓ Branch 2 taken 52 times.
64 if ( m_fstransient->iteration == 1) {
111
112 12 m_SimplifiedFSI->initStructure();
113
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 12 times.
12 if ( FelisceParam::verbose() > 0 )
114 PetscPrintf(MpiInfo::petscComm()," Structure initialized \n");
115
2/2
✓ Branch 0 taken 4 times.
✓ Branch 1 taken 8 times.
12 if (m_lumpedModelBC) {
116 // m_labelFluxListLumpedModelBC computation with LumpedModelBC labels
117 4 m_labelFluxListLumpedModelBC.clear();
118
2/2
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 4 times.
12 for(std::size_t i=0; i<m_lumpedModelBC->size(); i++) {
119 8 m_labelFluxListLumpedModelBC.push_back(FelisceParam::instance().lumpedModelBCLabel[i]);
120 }
121 }
122
123 12 m_bdfNS.initialize(m_linearProblem[m_idNSSimpl]->solution());
124 12 m_bdfNS.extrapolate(m_solExtrapolate);
125 12 m_linearProblem[m_idNSSimpl]->initExtrapol(m_solExtrapolate);
126 } else {
127 52 m_bdfNS.update(m_linearProblem[m_idNSSimpl]->solution());
128 52 m_bdfNS.extrapolate(m_solExtrapolate);
129 52 m_linearProblem[m_idNSSimpl]->updateExtrapol(m_solExtrapolate);
130 }
131
132 // lumpedModelBC bc computation
133
2/2
✓ Branch 0 taken 44 times.
✓ Branch 1 taken 20 times.
64 if (m_lumpedModelBC) {
134 44 m_SimplifiedFSI->handleWindkessel(m_fstransient->time);
135 // update the fluxes
136 44 m_linearProblem[m_idNSSimpl]->computeMeanQuantity(velocity,m_labelFluxListLumpedModelBC,m_fluxLumpedModelBC );
137
138
2/2
✓ Branch 1 taken 88 times.
✓ Branch 2 taken 44 times.
132 for(std::size_t i=0; i<m_lumpedModelBC->size(); i++) {
139 88 m_lumpedModelBC->Q(i) = m_fluxLumpedModelBC[i] ;
140 }
141 // update compliance if non linear compliance
142
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 44 times.
44 if ( FelisceParam::instance().lumpedModelBC_C_type == 2 ) {
143 static_cast<LinearProblemNSSimplifiedFSI*>(m_linearProblem[m_idNSSimpl])->updateVolume();
144 }
145 // iterate the lumpedModelBC parameters
146 44 m_lumpedModelBC->iterate();
147 }
148
149 64 m_bdfNS.computeRHSTime(m_fstransient->timeStep);
150 64 m_linearProblem[m_idNSSimpl]->gatherVectorBeforeAssembleMatrixRHS();
151 64 }
152
153 204 void NSSimplifiedFSIModel::solveNSSimplifiedFSI(FlagMatrixRHS flagMatrixRHS){
154
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 204 times.
204 if (FelisceParam::verbose() )
155 PetscPrintf(MpiInfo::petscComm(),"NS SimplifiedFSI \n");
156
157 204 m_linearProblem[m_idNSSimpl]->assembleMatrixRHS(MpiInfo::rankProc(),flagMatrixRHS); // assembly loop over the domain elements
158
3/4
✓ Branch 0 taken 140 times.
✓ Branch 1 taken 64 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 140 times.
204 if (flagMatrixRHS == FlagMatrixRHS::matrix_and_rhs || flagMatrixRHS == FlagMatrixRHS::only_matrix)
159 64 postAssemblingMatrixRHS(m_idNSSimpl); // sum "constant" matrix with time dependent
160
161
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 204 times.
204 if ( FelisceParam::verbose()>2)
162 std::cout<<"Time for assembling matrix and rhs: "<<chronoAssemble.diff()<<std::endl;
163
164
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 204 times.
204 if (FelisceParam::verbose()>2)
165 PetscPrintf(MpiInfo::petscComm(),"Applying BC \n");
166
167 204 m_linearProblem[m_idNSSimpl]->finalizeEssBCTransient();
168 204 m_linearProblem[m_idNSSimpl]->applyBC(FelisceParam::instance().essentialBoundaryConditionsMethod, MpiInfo::rankProc(), flagMatrixRHS,flagMatrixRHS);
169
170
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 204 times.
204 if (FelisceParam::verbose()>2)
171 PetscPrintf(MpiInfo::petscComm(),"Solving NS system \n");
172
173 204 m_linearProblem[m_idNSSimpl]->solve(MpiInfo::rankProc(), MpiInfo::numProc());
174 204 }
175
176 64 void NSSimplifiedFSIModel::finalizeForward() {
177 // update the structure! <with this kind of model, the structure is only a variable defined on the surface>
178
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 64 times.
64 if( FelisceParam::verbose() > 2 )
179 std::cout<<"["<<MpiInfo::rankProc()<<"] start - update structure"<<std::endl;
180 64 m_SimplifiedFSI->updateStructure();
181
182
1/2
✓ Branch 0 taken 64 times.
✗ Branch 1 not taken.
64 if( m_exportBdQuantities ) {
183 64 exportBdQuantities();
184 }
185
2/2
✓ Branch 1 taken 44 times.
✓ Branch 2 taken 20 times.
64 if( m_SimplifiedFSI->autoregulated())
186 44 m_SimplifiedFSI->updateControl(m_fstransient->time, MpiInfo::rankProc());
187 64 }
188
189 52 void NSSimplifiedFSIModel::forward() {
190
191 52 prepareForward();
192
193 52 solveNSSimplifiedFSI(FlagMatrixRHS::matrix_and_rhs);
194
195 52 finalizeForward();
196 52 }
197
198 64 void NSSimplifiedFSIModel::exportBdQuantities() {
199
200
3/6
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 64 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 64 times.
✗ Branch 9 not taken.
128 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"displacement"),
201 64 MpiInfo::rankProc(),
202 64 m_ios,
203
1/2
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
64 m_fstransient->time,
204 64 m_fstransient->iteration,
205
1/2
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
128 std::string("displacement"));
206
207
3/6
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 64 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 64 times.
✗ Branch 9 not taken.
128 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"normalField"),
208 64 MpiInfo::rankProc(),
209 64 m_ios,
210
1/2
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
64 m_fstransient->time,
211 64 m_fstransient->iteration,
212
1/2
✓ Branch 2 taken 64 times.
✗ Branch 3 not taken.
128 std::string("normal"));
213
1/2
✗ Branch 1 not taken.
✓ Branch 2 taken 64 times.
64 if ( FelisceParam::instance().exportDofPartition )
214 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"partition"),
215 MpiInfo::rankProc(),
216 m_ios,
217 m_fstransient->time,
218 m_fstransient->iteration,
219 std::string("partition"));
220
2/2
✓ Branch 1 taken 52 times.
✓ Branch 2 taken 12 times.
64 if (FelisceParam::instance().exportPrincipalDirectionsAndCurvature) {
221
3/6
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 52 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 52 times.
✗ Branch 9 not taken.
104 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::parallel,"meanCurvature"),
222 52 MpiInfo::rankProc(),
223 52 m_ios,
224
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
52 m_fstransient->time,
225 52 m_fstransient->iteration,
226
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
104 std::string("meanCurvature"));
227
228
3/6
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 52 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 52 times.
✗ Branch 9 not taken.
104 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"koiterCoefficients"),
229 52 MpiInfo::rankProc(),
230 52 m_ios,
231
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
52 m_fstransient->time,
232 52 m_fstransient->iteration,
233
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
104 std::string("koiterCoeff"));
234
235
3/6
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 52 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 52 times.
✗ Branch 9 not taken.
104 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"maxEigVec"),
236 52 MpiInfo::rankProc(),
237 52 m_ios,
238
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
52 m_fstransient->time,
239 52 m_fstransient->iteration,
240
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
104 std::string("maxEigenVector"));
241
242
3/6
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 52 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 52 times.
✗ Branch 9 not taken.
104 m_linearProblem[m_idNSSimpl]->writeSolutionFromVec(m_SimplifiedFSI->get(LinearProblemNSSimplifiedFSI::sequential,"minEigVec"),
243 52 MpiInfo::rankProc(),
244 52 m_ios,
245
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
52 m_fstransient->time,
246 52 m_fstransient->iteration,
247
1/2
✓ Branch 2 taken 52 times.
✗ Branch 3 not taken.
104 std::string("minEigenVector"));
248 }
249
250 64 }
251
252 }
253