Directory: | ./ |
---|---|
File: | Model/elasticityModel.cpp |
Date: | 2024-04-14 07:32:34 |
Exec | Total | Coverage | |
---|---|---|---|
Lines: | 45 | 55 | 81.8% |
Branches: | 7 | 20 | 35.0% |
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: S. Gilles | ||
13 | // | ||
14 | |||
15 | /*! | ||
16 | \file ElasticityModel.cpp | ||
17 | \authors S. Gilles | ||
18 | \date 21/02/2013 | ||
19 | \brief Model to manage elasticity problem. | ||
20 | */ | ||
21 | |||
22 | // System includes | ||
23 | |||
24 | // External includes | ||
25 | |||
26 | // Project includes | ||
27 | #include "Model/elasticityModel.hpp" | ||
28 | #include "Solver/linearProblemElasticity.hpp" | ||
29 | |||
30 | namespace felisce | ||
31 | { | ||
32 | |||
33 | 36 | ElasticityModel::ElasticityModel() | |
34 | 36 | : Model() | |
35 | { | ||
36 |
1/2✓ Branch 1 taken 36 times.
✗ Branch 2 not taken.
|
36 | m_name = "Elasticity"; |
37 | 36 | } | |
38 | |||
39 | /***********************************************************************************/ | ||
40 | /***********************************************************************************/ | ||
41 | |||
42 | 72 | ElasticityModel::~ElasticityModel() | |
43 | = default; | ||
44 | |||
45 | /***********************************************************************************/ | ||
46 | /***********************************************************************************/ | ||
47 | |||
48 | 36 | void ElasticityModel::SolveStaticProblem() | |
49 | { | ||
50 | #ifndef NDEBUG | ||
51 | static bool first_call = true; | ||
52 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 36 times.
|
36 | assert(first_call); |
53 | 36 | first_call = false; | |
54 | #endif // NDEBUG | ||
55 | |||
56 | 36 | PetscPrintf(MpiInfo::petscComm(), "\n----------------------------------------------\n"); | |
57 | 36 | PetscPrintf(MpiInfo::petscComm(), "Static problem\n"); | |
58 | 36 | PetscPrintf(MpiInfo::petscComm(), "----------------------------------------------\n"); | |
59 | |||
60 |
1/2✓ Branch 1 taken 36 times.
✗ Branch 2 not taken.
|
36 | LinearProblemElasticityDynamic& linear_problem = dynamic_cast<LinearProblemElasticityDynamic&>(*m_linearProblem[0]); |
61 | |||
62 | //Assembly loop on elements. | ||
63 | 36 | linear_problem.assembleMatrixRHS(MpiInfo::rankProc()); | |
64 | |||
65 | //Apply boundary conditions. | ||
66 | 36 | linear_problem.finalizeEssBCTransient(); | |
67 | 36 | linear_problem.applyBC(FelisceParam::instance().essentialBoundaryConditionsMethod, MpiInfo::rankProc(), FlagMatrixRHS::matrix_and_rhs, FlagMatrixRHS::matrix_and_rhs); | |
68 | |||
69 | // Solve the system. | ||
70 | 36 | linear_problem.solve(MpiInfo::rankProc(), MpiInfo::numProc()); | |
71 | |||
72 | // Update the displacement for the first time iteration. | ||
73 | 36 | linear_problem.UpdateCurrentDisplacement(); | |
74 | 36 | } | |
75 | |||
76 | /***********************************************************************************/ | ||
77 | /***********************************************************************************/ | ||
78 | |||
79 | 4 | void ElasticityModel::SolveDynamicProblem() | |
80 | { | ||
81 | 4 | SolveStaticProblem(); | |
82 | 4 | AssembleNewmarkMatrices(); | |
83 | |||
84 |
2/2✓ Branch 1 taken 120 times.
✓ Branch 2 taken 4 times.
|
124 | while (!hasFinished()) |
85 | 120 | forward(); | |
86 | 4 | } | |
87 | |||
88 | /***********************************************************************************/ | ||
89 | /***********************************************************************************/ | ||
90 | |||
91 | 120 | void ElasticityModel::updateTime(const FlagMatrixRHS flagMatrixRHS) | |
92 | { | ||
93 | IGNORE_UNUSED_FLAG_MATRIX_RHS; | ||
94 | 120 | m_fstransient->iteration++; | |
95 | 120 | m_fstransient->time +=m_fstransient->timeStep; | |
96 | 120 | } | |
97 | |||
98 | /***********************************************************************************/ | ||
99 | /***********************************************************************************/ | ||
100 | |||
101 | 4 | void ElasticityModel::AssembleNewmarkMatrices() | |
102 | { | ||
103 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | LinearProblemElasticityDynamic& linear_problem = dynamic_cast<LinearProblemElasticityDynamic&>(*m_linearProblem[0]); |
104 | 4 | linear_problem.GoToDynamic(); | |
105 | 4 | linear_problem.clearMatrixRHS(); | |
106 | 4 | linear_problem.assembleMatrixRHS(MpiInfo::rankProc(), FlagMatrixRHS::only_matrix); | |
107 | 4 | } | |
108 | |||
109 | /***********************************************************************************/ | ||
110 | /***********************************************************************************/ | ||
111 | |||
112 | 120 | void ElasticityModel::forward() | |
113 | { | ||
114 |
1/2✓ Branch 1 taken 120 times.
✗ Branch 2 not taken.
|
120 | LinearProblemElasticityDynamic& linear_problem = dynamic_cast<LinearProblemElasticityDynamic&>(*m_linearProblem[0]); |
115 | |||
116 | // Write solution for postprocessing (if required) | ||
117 | 120 | writeSolution(); | |
118 | |||
119 | // Advance time step. | ||
120 | 120 | updateTime(); | |
121 | |||
122 | // Print time information | ||
123 | 120 | printNewTimeIterationBanner(); | |
124 | |||
125 | // In the dynamic system only the RHS is modified each iteration; it is a mere algebric operation. | ||
126 | 120 | linear_problem.ComputeRHS(); | |
127 | |||
128 | //Apply boundary conditions. | ||
129 | 120 | linear_problem.finalizeEssBCTransient(); | |
130 | 120 | linear_problem.applyBC(FelisceParam::instance().essentialBoundaryConditionsMethod, MpiInfo::rankProc(), FlagMatrixRHS::matrix_and_rhs, | |
131 | FlagMatrixRHS::matrix_and_rhs, 0, false, ApplyNaturalBoundaryConditions::no); | ||
132 | |||
133 | //Solve linear system. | ||
134 | 120 | linear_problem.solve(MpiInfo::rankProc(), MpiInfo::numProc()); | |
135 | |||
136 | 120 | linear_problem.endIteration(); | |
137 | 120 | } | |
138 | |||
139 | /***********************************************************************************/ | ||
140 | /***********************************************************************************/ | ||
141 | |||
142 | ✗ | int ElasticityModel::getNstate() const | |
143 | { | ||
144 | ✗ | return m_linearProblem[0]->numDof(); | |
145 | } | ||
146 | |||
147 | /***********************************************************************************/ | ||
148 | /***********************************************************************************/ | ||
149 | |||
150 | ✗ | void ElasticityModel::getState(double* & state) | |
151 | { | ||
152 | ✗ | m_linearProblem[0]->getSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); | |
153 | } | ||
154 | |||
155 | /***********************************************************************************/ | ||
156 | /***********************************************************************************/ | ||
157 | |||
158 | ✗ | void ElasticityModel::getState_swig(double* data, felInt numDof) | |
159 | { | ||
160 | double * state; | ||
161 | ✗ | m_linearProblem[0]->getSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); | |
162 | |||
163 | ✗ | for (felInt i = 0; i < numDof ; i++) | |
164 | ✗ | data[i] = state[i]; | |
165 | } | ||
166 | |||
167 | /***********************************************************************************/ | ||
168 | /***********************************************************************************/ | ||
169 | |||
170 | ✗ | void ElasticityModel::setState(double* & state) | |
171 | { | ||
172 | ✗ | m_linearProblem[0]->setSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); | |
173 | } | ||
174 | |||
175 | |||
176 | } // namespace Elasticity | ||
177 |