suanPan
Loading...
Searching...
No Matches
Factory.hpp
Go to the documentation of this file.
1/*******************************************************************************
2 * Copyright (C) 2017-2023 Theodore Chang
3 *
4 * This program is free software: you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation, either version 3 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with this program. If not, see <http://www.gnu.org/licenses/>.
16 ******************************************************************************/
31#ifndef FACTORY_HPP
32#define FACTORY_HPP
33
34#include <future>
35#include <Toolbox/container.h>
36#include <Element/MappingDOF.h>
37#include <Domain/MetaMat/MetaMat>
38
39#ifdef SUANPAN_MAGMA
40#include <magmasparse.h>
41#endif
42
43enum class AnalysisType {
44 NONE,
45 DISP,
46 EIGEN,
47 BUCKLE,
48 STATICS,
50};
51
52enum class StorageScheme {
53 FULL,
54 BAND,
57 SPARSE,
59};
60
61enum class SolverType {
62 LAPACK,
63 SPIKE,
64 SUPERLU,
65 MUMPS,
66 CUDA,
67 PARDISO,
68 FGMRES,
69 MAGMA,
70 LIS
71};
72
73template<sp_d T> class Factory final {
74 unsigned n_size = 0; // number of degrees of freedom
75 unsigned n_lobw = 0; // low bandwidth
76 unsigned n_upbw = 0; // up bandwidth
77 unsigned n_sfbw = n_lobw + n_upbw; // matrix storage offset
78 unsigned n_rfld = 0; // reference load size
79 unsigned n_mpc = 0; // multipoint constraint size
80 uword n_elem = 0;
81
82 AnalysisType analysis_type = AnalysisType::NONE; // type of analysis
83 StorageScheme storage_type = StorageScheme::FULL; // type of analysis
84
85#ifdef SUANPAN_MAGMA
86 magma_dopts magma_setting{};
87#endif
88
89 bool nlgeom = false;
90
92 SolverSetting<T> setting{};
93
94 T error = 0.; // error produced by certain solvers
95
96 Col<T> ninja; // the result from A*X=B
97 Col<T> sushi; // modified right-hand side B
98
99 suanpan::set<uword> reference_dof;
100 SpMat<T> reference_load;
101
102 uvec auxiliary_encoding; // for constraints using multiplier method
103 Col<T> auxiliary_lambda; // for constraints using multiplier method
104 Col<T> auxiliary_resistance; // for constraints using multiplier method
105 Col<T> auxiliary_load; // for constraints using multiplier method
106 SpMat<T> auxiliary_stiffness; // for constraints using multiplier method
107
108 SpCol<T> trial_constraint_resistance;
109 SpCol<T> current_constraint_resistance;
110
111 T trial_time = T(0); // global trial (pseudo) time
112 T incre_time = T(0); // global incremental (pseudo) time
113 T current_time = T(0); // global current (pseudo) time
114 T pre_time = T(0); // global previous (pseudo) time
115
116 T strain_energy = T(0);
117 T kinetic_energy = T(0);
118 T viscous_energy = T(0);
119 T complementary_energy = T(0);
120 Col<T> momentum;
121
122 Col<T> trial_load_factor; // global trial load factor
123 Col<T> trial_load; // global trial load vector
124 Col<T> trial_settlement; // global trial displacement load vector
125 Col<T> trial_resistance; // global trial resistance vector
126 Col<T> trial_damping_force; // global trial damping force vector
127 Col<T> trial_inertial_force; // global trial inertial force vector
128 Col<T> trial_displacement; // global trial displacement vector
129 Col<T> trial_velocity; // global trial velocity vector
130 Col<T> trial_acceleration; // global trial acceleration vector
131 Col<T> trial_temperature; // global trial temperature vector
132
133 Col<T> incre_load_factor; // global incremental load vector
134 Col<T> incre_load; // global incremental load vector
135 Col<T> incre_settlement; // global incremental displacement load vector
136 Col<T> incre_resistance; // global incremental resistance vector
137 Col<T> incre_damping_force; // global incremental damping force vector
138 Col<T> incre_inertial_force; // global incremental inertial force vector
139 Col<T> incre_displacement; // global incremental displacement vector
140 Col<T> incre_velocity; // global incremental velocity vector
141 Col<T> incre_acceleration; // global incremental acceleration vector
142 Col<T> incre_temperature; // global incremental temperature vector
143
144 Col<T> current_load_factor; // global current load vector
145 Col<T> current_load; // global current load vector
146 Col<T> current_settlement; // global current displacement load vector
147 Col<T> current_resistance; // global current resistance vector
148 Col<T> current_damping_force; // global current damping force vector
149 Col<T> current_inertial_force; // global current inertial force vector
150 Col<T> current_displacement; // global current displacement vector
151 Col<T> current_velocity; // global current velocity vector
152 Col<T> current_acceleration; // global current acceleration vector
153 Col<T> current_temperature; // global current temperature vector
154
155 Col<T> pre_load_factor; // global previous load vector
156 Col<T> pre_load; // global previous load vector
157 Col<T> pre_settlement; // global previous displacement load vector
158 Col<T> pre_resistance; // global previous resistance vector
159 Col<T> pre_damping_force; // global previous damping force vector
160 Col<T> pre_inertial_force; // global previous inertial force vector
161 Col<T> pre_displacement; // global previous displacement vector
162 Col<T> pre_velocity; // global previous velocity vector
163 Col<T> pre_acceleration; // global previous acceleration vector
164 Col<T> pre_temperature; // global previous temperature vector
165
166 shared_ptr<MetaMat<T>> global_mass = nullptr; // global mass matrix
167 shared_ptr<MetaMat<T>> global_damping = nullptr; // global damping matrix
168 shared_ptr<MetaMat<T>> global_stiffness = nullptr; // global stiffness matrix
169 shared_ptr<MetaMat<T>> global_geometry = nullptr; // global geometry matrix
170
171 std::vector<std::mutex> global_mutex = std::vector<std::mutex>(20);
172
173 Col<T> eigenvalue; // eigenvalues
174
175 Mat<T> eigenvector; // eigenvectors
176
177 unique_ptr<MetaMat<T>> get_basic_container();
178 unique_ptr<MetaMat<T>> get_matrix_container();
179
180 void assemble_matrix_helper(shared_ptr<MetaMat<T>>&, const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
181
182public:
183 const bool initialized = false;
184
186
187 void set_size(unsigned);
188 [[nodiscard]] unsigned get_size() const;
189
190 void set_entry(uword);
191 [[nodiscard]] uword get_entry() const;
192
193 void set_nlgeom(bool);
194 [[nodiscard]] bool is_nlgeom() const;
195
197 [[nodiscard]] SolverType get_solver_type() const;
198
200 [[nodiscard]] const SolverSetting<double>& get_solver_setting() const;
201
202#ifdef SUANPAN_MAGMA
203 void set_solver_setting(const magma_dopts& magma_opt) { magma_setting = magma_opt; }
204
205 [[nodiscard]] const magma_dopts& get_magma_setting() const { return magma_setting; }
206#endif
207
209 [[nodiscard]] AnalysisType get_analysis_type() const;
210
212 [[nodiscard]] StorageScheme get_storage_scheme() const;
213
214 [[nodiscard]] bool is_sparse() const;
215
216 void set_bandwidth(unsigned, unsigned);
217 void get_bandwidth(unsigned&, unsigned&) const;
218
220 void set_reference_size(unsigned);
221 [[nodiscard]] unsigned get_reference_size() const;
222
223 void update_reference_dof(const uvec&);
225 [[nodiscard]] const suanpan::set<uword>& get_reference_dof() const;
226
227 void set_error(T);
228 T get_error() const;
229
230 /*************************INITIALIZER*************************/
231
232 int initialize();
233
235 void initialize_load();
241 void initialize_velocity();
245
246 void initialize_mass();
247 void initialize_damping();
249 void initialize_geometry();
250 void initialize_eigen();
251
252 /*************************SETTER*************************/
253
254 void set_ninja(const Col<T>&);
255 void set_sushi(const Col<T>&);
256
257 void update_sushi_by(const Col<T>&);
258
259 void set_mpc(unsigned);
260
261 void set_reference_load(const SpMat<T>&);
262
263 void set_trial_time(T);
264 void set_trial_load_factor(const Col<T>&);
265 void set_trial_load(const Col<T>&);
266 void set_trial_settlement(const Col<T>&);
267 void set_trial_resistance(const Col<T>&);
268 void set_trial_damping_force(const Col<T>&);
269 void set_trial_inertial_force(const Col<T>&);
270 void set_trial_displacement(const Col<T>&);
271 void set_trial_velocity(const Col<T>&);
272 void set_trial_acceleration(const Col<T>&);
273 void set_trial_temperature(const Col<T>&);
274
275 void set_incre_time(T);
276 void set_incre_load_factor(const Col<T>&);
277 void set_incre_load(const Col<T>&);
278 void set_incre_settlement(const Col<T>&);
279 void set_incre_resistance(const Col<T>&);
280 void set_incre_damping_force(const Col<T>&);
281 void set_incre_inertial_force(const Col<T>&);
282 void set_incre_displacement(const Col<T>&);
283 void set_incre_velocity(const Col<T>&);
284 void set_incre_acceleration(const Col<T>&);
285 void set_incre_temperature(const Col<T>&);
286
287 void set_current_time(T);
288 void set_current_load_factor(const Col<T>&);
289 void set_current_load(const Col<T>&);
290 void set_current_settlement(const Col<T>&);
291 void set_current_resistance(const Col<T>&);
292 void set_current_damping_force(const Col<T>&);
293 void set_current_inertial_force(const Col<T>&);
294 void set_current_displacement(const Col<T>&);
295 void set_current_velocity(const Col<T>&);
296 void set_current_acceleration(const Col<T>&);
297 void set_current_temperature(const Col<T>&);
298
299 void set_pre_time(T);
300 void set_pre_load_factor(const Col<T>&);
301 void set_pre_load(const Col<T>&);
302 void set_pre_settlement(const Col<T>&);
303 void set_pre_resistance(const Col<T>&);
304 void set_pre_damping_force(const Col<T>&);
305 void set_pre_inertial_force(const Col<T>&);
306 void set_pre_displacement(const Col<T>&);
307 void set_pre_velocity(const Col<T>&);
308 void set_pre_acceleration(const Col<T>&);
309 void set_pre_temperature(const Col<T>&);
310
311 void set_mass(const shared_ptr<MetaMat<T>>&);
312 void set_damping(const shared_ptr<MetaMat<T>>&);
313 void set_stiffness(const shared_ptr<MetaMat<T>>&);
314 void set_geometry(const shared_ptr<MetaMat<T>>&);
315
316 void set_eigenvalue(const Col<T>&);
317 void set_eigenvector(const Mat<T>&);
318
319 /*************************GETTER*************************/
320
321 const Col<T>& get_ninja() const;
322 const Col<T>& get_sushi() const;
323
324 [[nodiscard]] unsigned get_mpc() const;
325
326 const SpMat<T>& get_reference_load() const;
327
328 [[nodiscard]] const uvec& get_auxiliary_encoding() const;
329 const Col<T>& get_auxiliary_lambda() const;
330 const Col<T>& get_auxiliary_resistance() const;
331 const Col<T>& get_auxiliary_load() const;
332 const SpMat<T>& get_auxiliary_stiffness() const;
333
334 const SpCol<T>& get_trial_constraint_resistance() const;
335 const SpCol<T>& get_current_constraint_resistance() const;
336
341 const Col<T>& get_momentum();
342
343 T get_trial_time() const;
344 const Col<T>& get_trial_load_factor() const;
345 const Col<T>& get_trial_load() const;
346 const Col<T>& get_trial_settlement() const;
347 const Col<T>& get_trial_resistance() const;
348 const Col<T>& get_trial_damping_force() const;
349 const Col<T>& get_trial_inertial_force() const;
350 const Col<T>& get_trial_displacement() const;
351 const Col<T>& get_trial_velocity() const;
352 const Col<T>& get_trial_acceleration() const;
353 const Col<T>& get_trial_temperature() const;
354
355 T get_incre_time() const;
356 const Col<T>& get_incre_load_factor() const;
357 const Col<T>& get_incre_load() const;
358 const Col<T>& get_incre_settlement() const;
359 const Col<T>& get_incre_resistance() const;
360 const Col<T>& get_incre_damping_force() const;
361 const Col<T>& get_incre_inertial_force() const;
362 const Col<T>& get_incre_displacement() const;
363 const Col<T>& get_incre_velocity() const;
364 const Col<T>& get_incre_acceleration() const;
365 const Col<T>& get_incre_temperature() const;
366
367 T get_current_time() const;
368 const Col<T>& get_current_load_factor() const;
369 const Col<T>& get_current_load() const;
370 const Col<T>& get_current_settlement() const;
371 const Col<T>& get_current_resistance() const;
372 const Col<T>& get_current_damping_force() const;
373 const Col<T>& get_current_inertial_force() const;
374 const Col<T>& get_current_displacement() const;
375 const Col<T>& get_current_velocity() const;
376 const Col<T>& get_current_acceleration() const;
377 const Col<T>& get_current_temperature() const;
378
379 T get_pre_time() const;
380 const Col<T>& get_pre_load_factor() const;
381 const Col<T>& get_pre_load() const;
382 const Col<T>& get_pre_settlement() const;
383 const Col<T>& get_pre_resistance() const;
384 const Col<T>& get_pre_damping_force() const;
385 const Col<T>& get_pre_inertial_force() const;
386 const Col<T>& get_pre_displacement() const;
387 const Col<T>& get_pre_velocity() const;
388 const Col<T>& get_pre_acceleration() const;
389 const Col<T>& get_pre_temperature() const;
390
391 const shared_ptr<MetaMat<T>>& get_mass() const;
392 const shared_ptr<MetaMat<T>>& get_damping() const;
393 const shared_ptr<MetaMat<T>>& get_stiffness() const;
394 const shared_ptr<MetaMat<T>>& get_geometry() const;
395
396 std::mutex& get_auxiliary_encoding_mutex();
397 std::mutex& get_auxiliary_resistance_mutex();
398 std::mutex& get_auxiliary_load_mutex();
399 std::mutex& get_auxiliary_stiffness_mutex();
400
402
403 std::mutex& get_trial_load_mutex();
404 std::mutex& get_trial_settlement_mutex();
405
406 std::mutex& get_mass_mutex();
407 std::mutex& get_damping_mutex();
408 std::mutex& get_stiffness_mutex();
409 std::mutex& get_geometry_mutex();
410
411 const Col<T>& get_eigenvalue() const;
412 const Mat<T>& get_eigenvector() const;
413
414 /*************************UPDATER*************************/
415
416 void update_trial_time(T);
417 void update_trial_load_factor(const Col<T>&);
418 void update_trial_load(const Col<T>&);
419 void update_trial_settlement(const Col<T>&);
420 void update_trial_resistance(const Col<T>&);
421 void update_trial_damping_force(const Col<T>&);
422 void update_trial_inertial_force(const Col<T>&);
423 void update_trial_displacement(const Col<T>&);
424 void update_trial_velocity(const Col<T>&);
425 void update_trial_acceleration(const Col<T>&);
426 void update_trial_temperature(const Col<T>&);
427
428 void update_incre_time(T);
429 void update_incre_load_factor(const Col<T>&);
430 void update_incre_load(const Col<T>&);
431 void update_incre_settlement(const Col<T>&);
432 void update_incre_resistance(const Col<T>&);
433 void update_incre_damping_force(const Col<T>&);
434 void update_incre_inertial_force(const Col<T>&);
435 void update_incre_displacement(const Col<T>&);
436 void update_incre_velocity(const Col<T>&);
437 void update_incre_acceleration(const Col<T>&);
438 void update_incre_temperature(const Col<T>&);
439
441 void update_current_load_factor(const Col<T>&);
442 void update_current_load(const Col<T>&);
443 void update_current_settlement(const Col<T>&);
444 void update_current_resistance(const Col<T>&);
445 void update_current_damping_force(const Col<T>&);
446 void update_current_inertial_force(const Col<T>&);
447 void update_current_displacement(const Col<T>&);
448 void update_current_velocity(const Col<T>&);
449 void update_current_acceleration(const Col<T>&);
450 void update_current_temperature(const Col<T>&);
451
453 void update_trial_load_factor_by(const Col<T>&);
454 void update_trial_load_by(const Col<T>&);
455 void update_trial_settlement_by(const Col<T>&);
456 void update_trial_resistance_by(const Col<T>&);
457 void update_trial_damping_force_by(const Col<T>&);
458 void update_trial_inertial_force_by(const Col<T>&);
459 void update_trial_displacement_by(const Col<T>&);
460 void update_trial_velocity_by(const Col<T>&);
461 void update_trial_acceleration_by(const Col<T>&);
462 void update_trial_temperature_by(const Col<T>&);
463
465 void update_incre_load_factor_by(const Col<T>&);
466 void update_incre_load_by(const Col<T>&);
467 void update_incre_settlement_by(const Col<T>&);
468 void update_incre_resistance_by(const Col<T>&);
469 void update_incre_damping_force_by(const Col<T>&);
470 void update_incre_inertial_force_by(const Col<T>&);
471 void update_incre_displacement_by(const Col<T>&);
472 void update_incre_velocity_by(const Col<T>&);
473 void update_incre_acceleration_by(const Col<T>&);
474 void update_incre_temperature_by(const Col<T>&);
475
477 void update_current_load_factor_by(const Col<T>&);
478 void update_current_load_by(const Col<T>&);
479 void update_current_settlement_by(const Col<T>&);
480 void update_current_resistance_by(const Col<T>&);
481 void update_current_damping_force_by(const Col<T>&);
482 void update_current_inertial_force_by(const Col<T>&);
483 void update_current_displacement_by(const Col<T>&);
484 void update_current_velocity_by(const Col<T>&);
485 void update_current_acceleration_by(const Col<T>&);
486 void update_current_temperature_by(const Col<T>&);
487
488 /*************************FRIEND*************************/
489
490 Col<T>& modify_ninja();
491 Col<T>& modify_sushi();
492
494 SpMat<T>& modify_reference_load();
495
497 Col<T>& modify_auxiliary_lambda();
499 Col<T>& modify_auxiliary_load();
500 SpMat<T>& modify_auxiliary_stiffness();
501
504
506 Col<T>& modify_trial_load_factor();
507 Col<T>& modify_trial_load();
508 Col<T>& modify_trial_settlement();
509 Col<T>& modify_trial_resistance();
513 Col<T>& modify_trial_velocity();
515 Col<T>& modify_trial_temperature();
516
518 Col<T>& modify_incre_load_factor();
519 Col<T>& modify_incre_load();
520 Col<T>& modify_incre_settlement();
521 Col<T>& modify_incre_resistance();
525 Col<T>& modify_incre_velocity();
527 Col<T>& modify_incre_temperature();
528
531 Col<T>& modify_current_load();
537 Col<T>& modify_current_velocity();
540
542 Col<T>& modify_pre_load_factor();
543 Col<T>& modify_pre_load();
544 Col<T>& modify_pre_settlement();
545 Col<T>& modify_pre_resistance();
546 Col<T>& modify_pre_damping_force();
548 Col<T>& modify_pre_displacement();
549 Col<T>& modify_pre_velocity();
550 Col<T>& modify_pre_acceleration();
551 Col<T>& modify_pre_temperature();
552
553 shared_ptr<MetaMat<T>>& modify_mass();
554 shared_ptr<MetaMat<T>>& modify_damping();
555 shared_ptr<MetaMat<T>>& modify_stiffness();
556 shared_ptr<MetaMat<T>>& modify_geometry();
557
558 Col<T>& modify_eigenvalue();
559 Mat<T>& modify_eigenvector();
560
561 /*************************STATUS*************************/
562
563 void commit_energy();
564 void clear_energy();
565
566 void commit_status();
567 void commit_time();
568 void commit_load_factor();
569 void commit_load();
570 void commit_settlement();
571 void commit_resistance();
574 void commit_displacement();
575 void commit_velocity();
576 void commit_acceleration();
577 void commit_temperature();
579
580 void commit_pre_status();
581 void commit_pre_time();
583 void commit_pre_load();
589 void commit_pre_velocity();
592
593 void clear_status();
594 void clear_time();
595 void clear_load_factor();
596 void clear_load();
597 void clear_settlement();
598 void clear_resistance();
599 void clear_damping_force();
601 void clear_displacement();
602 void clear_velocity();
603 void clear_acceleration();
604 void clear_temperature();
606
607 void reset_status();
608 void reset_time();
609 void reset_load_factor();
610 void reset_load();
611 void reset_settlement();
612 void reset_resistance();
613 void reset_damping_force();
615 void reset_displacement();
616 void reset_velocity();
617 void reset_acceleration();
618 void reset_temperature();
620
621 void clear_eigen();
622 void clear_mass();
623 void clear_damping();
624 void clear_stiffness();
625 void clear_geometry();
626 void clear_auxiliary();
627
628 void reset();
629
630 /*************************ASSEMBLER*************************/
631
632 void assemble_resistance(const Mat<T>&, const uvec&);
633 void assemble_damping_force(const Mat<T>&, const uvec&);
634 void assemble_inertial_force(const Mat<T>&, const uvec&);
635
636 void assemble_mass(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
637 void assemble_damping(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
638 void assemble_stiffness(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
639 void assemble_geometry(const Mat<T>&, const uvec&, const std::vector<MappingDOF>&);
640
641 void assemble_stiffness(const SpMat<T>&, const uvec&);
642
643 /*************************UTILITY*************************/
644
645 void print() const;
646};
647
648template<sp_d T> Factory<T>::Factory(const unsigned D, const AnalysisType AT, const StorageScheme SS)
649 : n_size(D)
650 , analysis_type(AT)
651 , storage_type(SS) {}
652
653template<sp_d T> void Factory<T>::set_size(const unsigned D) {
654 if(D == n_size) return;
655 n_size = D;
656 access::rw(initialized) = false;
657}
658
659template<sp_d T> unsigned Factory<T>::get_size() const { return n_size; }
660
661template<sp_d T> void Factory<T>::set_entry(const uword N) {
662 n_elem = N;
663 if(n_elem > std::numeric_limits<int>::max()) throw invalid_argument("too many elements");
664}
665
666template<sp_d T> uword Factory<T>::get_entry() const { return n_elem; }
667
668template<sp_d T> void Factory<T>::set_nlgeom(const bool B) {
669 if(B == nlgeom) return;
670 nlgeom = B;
671 access::rw(initialized) = false;
672}
673
674template<sp_d T> bool Factory<T>::is_nlgeom() const { return nlgeom; }
675
676template<sp_d T> void Factory<T>::set_solver_type(const SolverType E) { solver = E; }
677
678template<sp_d T> SolverType Factory<T>::get_solver_type() const { return solver; }
679
680template<sp_d T> void Factory<T>::set_solver_setting(const SolverSetting<double>& SS) { setting = SS; }
681
682template<sp_d T> const SolverSetting<double>& Factory<T>::get_solver_setting() const { return setting; }
683
684template<sp_d T> void Factory<T>::set_analysis_type(const AnalysisType AT) {
685 if(AT == analysis_type) return;
686 analysis_type = AT;
687 access::rw(initialized) = false;
688}
689
690template<sp_d T> AnalysisType Factory<T>::get_analysis_type() const { return analysis_type; }
691
692template<sp_d T> void Factory<T>::set_storage_scheme(const StorageScheme SS) {
693 if(SS == storage_type) return;
694 storage_type = SS;
695 access::rw(initialized) = false;
696}
697
698template<sp_d T> StorageScheme Factory<T>::get_storage_scheme() const { return storage_type; }
699
700template<sp_d T> bool Factory<T>::is_sparse() const { return StorageScheme::SPARSE == storage_type || StorageScheme::SPARSESYMM == storage_type; }
701
702template<sp_d T> void Factory<T>::set_bandwidth(const unsigned L, const unsigned U) {
703 if(L == n_lobw && U == n_upbw) return;
704 n_lobw = L;
705 n_upbw = U;
706 n_sfbw = L + U;
707 access::rw(initialized) = false;
708}
709
710template<sp_d T> void Factory<T>::get_bandwidth(unsigned& L, unsigned& U) const {
711 L = n_lobw;
712 U = n_upbw;
713}
714
715template<sp_d T> void Factory<T>::update_reference_size() { n_rfld = static_cast<unsigned>(reference_dof.size()); }
716
717template<sp_d T> void Factory<T>::set_reference_size(const unsigned S) {
718 if(S == n_rfld) return;
719 n_rfld = S;
720}
721
722template<sp_d T> unsigned Factory<T>::get_reference_size() const { return n_rfld; }
723
724template<sp_d T> void Factory<T>::update_reference_dof(const uvec& S) { reference_dof.insert(S.cbegin(), S.cend()); }
725
726template<sp_d T> void Factory<T>::set_reference_dof(const suanpan::set<uword>& D) { reference_dof = D; }
727
728template<sp_d T> const suanpan::set<uword>& Factory<T>::get_reference_dof() const { return reference_dof; }
729
730template<sp_d T> void Factory<T>::set_error(const T E) { error = E; }
731
732template<sp_d T> T Factory<T>::get_error() const { return error; }
733
734template<sp_d T> int Factory<T>::initialize() {
735 reference_dof.clear(); // clear reference dof vector in every step
736
737 if(initialized || n_size == 0) return 0;
738
739 ninja.zeros(n_size);
740 sushi.zeros(n_size);
741
742 reset();
743
744 switch(analysis_type) {
746 initialize_displacement();
747 break;
749 initialize_mass();
750 initialize_stiffness();
751 initialize_eigen();
752 break;
755 initialize_load();
756 initialize_resistance();
757 initialize_displacement();
758 initialize_stiffness();
759 initialize_geometry();
760 break;
762 initialize_load();
763 initialize_resistance();
764 initialize_damping_force();
765 initialize_inertial_force();
766 initialize_displacement();
767 initialize_velocity();
768 initialize_acceleration();
769 initialize_mass();
770 initialize_damping();
771 initialize_stiffness();
772 initialize_geometry();
773 break;
775 break;
776 }
777
778 initialize_auxiliary_resistance();
779
780 access::rw(initialized) = true;
781
782 return 0;
783}
784
785template<sp_d T> void Factory<T>::initialize_load_factor() {
786 if(n_rfld == 0) return;
787
788 trial_load_factor.zeros(n_rfld);
789 incre_load_factor.zeros(n_rfld);
790 current_load_factor.zeros(n_rfld);
791
792 reference_load.zeros(n_size, n_rfld);
793}
794
795template<sp_d T> void Factory<T>::initialize_load() {
796 trial_load.zeros(n_size);
797 incre_load.zeros(n_size);
798 current_load.zeros(n_size);
799}
800
801template<sp_d T> void Factory<T>::initialize_settlement() {
802 trial_settlement.zeros(n_size);
803 incre_settlement.zeros(n_size);
804 current_settlement.zeros(n_size);
805}
806
807template<sp_d T> void Factory<T>::initialize_resistance() {
808 trial_resistance.zeros(n_size);
809 incre_resistance.zeros(n_size);
810 current_resistance.zeros(n_size);
811}
812
814 trial_damping_force.zeros(n_size);
815 incre_damping_force.zeros(n_size);
816 current_damping_force.zeros(n_size);
817}
818
820 trial_inertial_force.zeros(n_size);
821 incre_inertial_force.zeros(n_size);
822 current_inertial_force.zeros(n_size);
823}
824
825template<sp_d T> void Factory<T>::initialize_displacement() {
826 trial_displacement.zeros(n_size);
827 incre_displacement.zeros(n_size);
828 current_displacement.zeros(n_size);
829}
830
831template<sp_d T> void Factory<T>::initialize_velocity() {
832 trial_velocity.zeros(n_size);
833 incre_velocity.zeros(n_size);
834 current_velocity.zeros(n_size);
835}
836
837template<sp_d T> void Factory<T>::initialize_acceleration() {
838 trial_acceleration.zeros(n_size);
839 incre_acceleration.zeros(n_size);
840 current_acceleration.zeros(n_size);
841}
842
843template<sp_d T> void Factory<T>::initialize_temperature() {
844 trial_temperature.zeros(n_size);
845 incre_temperature.zeros(n_size);
846 current_temperature.zeros(n_size);
847}
848
850 trial_constraint_resistance.zeros(n_size);
851 current_constraint_resistance.zeros(n_size);
852}
853
854template<sp_d T> void Factory<T>::initialize_mass() { global_mass = get_matrix_container(); }
855
856template<sp_d T> void Factory<T>::initialize_damping() { global_damping = get_matrix_container(); }
857
858template<sp_d T> void Factory<T>::initialize_stiffness() { global_stiffness = get_matrix_container(); }
859
860template<sp_d T> void Factory<T>::initialize_geometry() {
861 if(!nlgeom) return;
862
863 global_geometry = get_matrix_container();
864}
865
866template<sp_d T> void Factory<T>::initialize_eigen() {
867 eigenvalue.zeros(n_size);
868 eigenvector.zeros(n_size, n_size);
869}
870
871template<sp_d T> void Factory<T>::set_ninja(const Col<T>& N) { ninja = N; }
872
873template<sp_d T> void Factory<T>::set_sushi(const Col<T>& S) { sushi = S; }
874
875template<sp_d T> void Factory<T>::update_sushi_by(const Col<T>& S) { sushi += S; }
876
877template<sp_d T> void Factory<T>::set_mpc(const unsigned S) {
878 n_mpc = S;
879 auxiliary_encoding.zeros(n_mpc);
880 auxiliary_resistance.zeros(n_mpc);
881 auxiliary_load.zeros(n_mpc);
882 auxiliary_stiffness.zeros(n_size, n_mpc);
883}
884
885template<sp_d T> void Factory<T>::set_reference_load(const SpMat<T>& L) { reference_load = L; }
886
887template<sp_d T> void Factory<T>::set_mass(const shared_ptr<MetaMat<T>>& M) { global_mass = M; }
888
889template<sp_d T> void Factory<T>::set_damping(const shared_ptr<MetaMat<T>>& C) { global_damping = C; }
890
891template<sp_d T> void Factory<T>::set_stiffness(const shared_ptr<MetaMat<T>>& K) { global_stiffness = K; }
892
893template<sp_d T> void Factory<T>::set_geometry(const shared_ptr<MetaMat<T>>& G) { global_geometry = G; }
894
895template<sp_d T> void Factory<T>::set_eigenvalue(const Col<T>& L) { eigenvalue = L; }
896
897template<sp_d T> void Factory<T>::set_eigenvector(const Mat<T>& V) { eigenvector = V; }
898
899template<sp_d T> const Col<T>& Factory<T>::get_ninja() const { return ninja; }
900
901template<sp_d T> const Col<T>& Factory<T>::get_sushi() const { return sushi; }
902
903template<sp_d T> unsigned Factory<T>::get_mpc() const { return n_mpc; }
904
905template<sp_d T> const SpMat<T>& Factory<T>::get_reference_load() const { return reference_load; }
906
907template<sp_d T> const uvec& Factory<T>::get_auxiliary_encoding() const { return auxiliary_encoding; }
908
909template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_lambda() const { return auxiliary_lambda; }
910
911template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_resistance() const { return auxiliary_resistance; }
912
913template<sp_d T> const Col<T>& Factory<T>::get_auxiliary_load() const { return auxiliary_load; }
914
915template<sp_d T> const SpMat<T>& Factory<T>::get_auxiliary_stiffness() const { return auxiliary_stiffness; }
916
917template<sp_d T> const SpCol<T>& Factory<T>::get_trial_constraint_resistance() const { return trial_constraint_resistance; }
918
919template<sp_d T> const SpCol<T>& Factory<T>::get_current_constraint_resistance() const { return current_constraint_resistance; }
920
921template<sp_d T> T Factory<T>::get_strain_energy() { return strain_energy; }
922
923template<sp_d T> T Factory<T>::get_kinetic_energy() { return kinetic_energy; }
924
925template<sp_d T> T Factory<T>::get_viscous_energy() { return viscous_energy; }
926
927template<sp_d T> T Factory<T>::get_complementary_energy() { return complementary_energy; }
928
929template<sp_d T> const Col<T>& Factory<T>::get_momentum() { return momentum; }
930
931template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_mass() const { return global_mass; }
932
933template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_damping() const { return global_damping; }
934
935template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_stiffness() const { return global_stiffness; }
936
937template<sp_d T> const shared_ptr<MetaMat<T>>& Factory<T>::get_geometry() const { return global_geometry; }
938
939template<sp_d T> std::mutex& Factory<T>::get_auxiliary_encoding_mutex() { return global_mutex.at(0); }
940
941template<sp_d T> std::mutex& Factory<T>::get_auxiliary_resistance_mutex() { return global_mutex.at(1); }
942
943template<sp_d T> std::mutex& Factory<T>::get_auxiliary_load_mutex() { return global_mutex.at(2); }
944
945template<sp_d T> std::mutex& Factory<T>::get_auxiliary_stiffness_mutex() { return global_mutex.at(3); }
946
947template<sp_d T> std::mutex& Factory<T>::get_trial_constraint_resistance_mutex() { return global_mutex.at(4); }
948
949template<sp_d T> std::mutex& Factory<T>::get_trial_load_mutex() { return global_mutex.at(5); }
950
951template<sp_d T> std::mutex& Factory<T>::get_trial_settlement_mutex() { return global_mutex.at(6); }
952
953template<sp_d T> std::mutex& Factory<T>::get_mass_mutex() { return global_mutex.at(7); }
954
955template<sp_d T> std::mutex& Factory<T>::get_damping_mutex() { return global_mutex.at(8); }
956
957template<sp_d T> std::mutex& Factory<T>::get_stiffness_mutex() { return global_mutex.at(9); }
958
959template<sp_d T> std::mutex& Factory<T>::get_geometry_mutex() { return global_mutex.at(10); }
960
961template<sp_d T> const Col<T>& Factory<T>::get_eigenvalue() const { return eigenvalue; }
962
963template<sp_d T> const Mat<T>& Factory<T>::get_eigenvector() const { return eigenvector; }
964
965template<sp_d T> void Factory<T>::commit_energy() {
966 auto se = std::async([&] { if(!trial_resistance.empty() && !incre_displacement.empty()) strain_energy += .5 * dot(trial_resistance + current_resistance, incre_displacement); });
967 auto ke = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) kinetic_energy = .5 * dot(global_mass * trial_velocity, trial_velocity); });
968 auto ve = std::async([&] { if(!trial_damping_force.empty() && !incre_displacement.empty()) viscous_energy += .5 * dot(trial_damping_force + current_damping_force, incre_displacement); });
969 auto ce = std::async([&] { if(!trial_displacement.empty() && !incre_resistance.empty()) complementary_energy += .5 * dot(trial_displacement + current_displacement, incre_resistance); });
970 auto mm = std::async([&] { if(!trial_inertial_force.empty() && !trial_velocity.empty()) momentum = global_mass * trial_velocity; });
971
972 se.get();
973 ke.get();
974 ve.get();
975 ce.get();
976 mm.get();
977}
978
979template<sp_d T> void Factory<T>::clear_energy() {
980 strain_energy = T(0);
981 kinetic_energy = T(0);
982 viscous_energy = T(0);
983 complementary_energy = T(0);
984 momentum.zeros();
985}
986
987template<sp_d T> void Factory<T>::commit_status() {
988 ninja.zeros();
989
990 commit_energy();
991
992 commit_time();
993 commit_load_factor();
994 commit_load();
995 commit_settlement();
996 commit_resistance();
997 commit_damping_force();
998 commit_inertial_force();
999 commit_displacement();
1000 commit_velocity();
1001 commit_acceleration();
1002 commit_temperature();
1003 commit_auxiliary_resistance();
1004}
1005
1006template<sp_d T> void Factory<T>::commit_time() {
1007 current_time = trial_time;
1008 incre_time = T(0);
1009}
1010
1011template<sp_d T> void Factory<T>::commit_load_factor() {
1012 if(trial_load_factor.is_empty()) return;
1013 current_load_factor = trial_load_factor;
1014 incre_load_factor.zeros();
1015}
1016
1017template<sp_d T> void Factory<T>::commit_load() {
1018 if(trial_load.is_empty()) return;
1019 current_load = trial_load;
1020 incre_load.zeros();
1021}
1022
1023template<sp_d T> void Factory<T>::commit_settlement() {
1024 if(trial_settlement.is_empty()) return;
1025 current_settlement = trial_settlement;
1026 incre_settlement.zeros();
1027}
1028
1029template<sp_d T> void Factory<T>::commit_resistance() {
1030 if(trial_resistance.is_empty()) return;
1031 current_resistance = trial_resistance;
1032 incre_resistance.zeros();
1033}
1034
1035template<sp_d T> void Factory<T>::commit_damping_force() {
1036 if(trial_damping_force.is_empty()) return;
1037 current_damping_force = trial_damping_force;
1038 incre_damping_force.zeros();
1039}
1040
1041template<sp_d T> void Factory<T>::commit_inertial_force() {
1042 if(trial_inertial_force.is_empty()) return;
1043 current_inertial_force = trial_inertial_force;
1044 incre_inertial_force.zeros();
1045}
1046
1047template<sp_d T> void Factory<T>::commit_displacement() {
1048 if(trial_displacement.is_empty()) return;
1049 current_displacement = trial_displacement;
1050 incre_displacement.zeros();
1051}
1052
1053template<sp_d T> void Factory<T>::commit_velocity() {
1054 if(trial_velocity.is_empty()) return;
1055 current_velocity = trial_velocity;
1056 incre_velocity.zeros();
1057}
1058
1059template<sp_d T> void Factory<T>::commit_acceleration() {
1060 if(trial_acceleration.is_empty()) return;
1061 current_acceleration = trial_acceleration;
1062 incre_acceleration.zeros();
1063}
1064
1065template<sp_d T> void Factory<T>::commit_temperature() {
1066 if(trial_temperature.is_empty()) return;
1067 current_temperature = trial_temperature;
1068 incre_temperature.zeros();
1069}
1070
1072 if(trial_constraint_resistance.is_empty()) return;
1073 current_constraint_resistance = trial_constraint_resistance;
1074}
1075
1076template<sp_d T> void Factory<T>::commit_pre_status() {
1077 commit_pre_time();
1078 commit_pre_load_factor();
1079 commit_pre_load();
1080 commit_pre_settlement();
1081 commit_pre_resistance();
1082 commit_pre_damping_force();
1083 commit_pre_inertial_force();
1084 commit_pre_displacement();
1085 commit_pre_velocity();
1086 commit_pre_acceleration();
1087 commit_pre_temperature();
1088}
1089
1090template<sp_d T> void Factory<T>::commit_pre_time() { pre_time = current_time; }
1091
1092template<sp_d T> void Factory<T>::commit_pre_load_factor() { if(!current_load_factor.is_empty()) pre_load_factor = current_load_factor; }
1093
1094template<sp_d T> void Factory<T>::commit_pre_load() { if(!current_load.is_empty()) pre_load = current_load; }
1095
1096template<sp_d T> void Factory<T>::commit_pre_settlement() { if(!current_settlement.is_empty()) pre_settlement = current_settlement; }
1097
1098template<sp_d T> void Factory<T>::commit_pre_resistance() { if(!current_resistance.is_empty()) pre_resistance = current_resistance; }
1099
1100template<sp_d T> void Factory<T>::commit_pre_damping_force() { if(!current_damping_force.is_empty()) pre_damping_force = current_damping_force; }
1101
1102template<sp_d T> void Factory<T>::commit_pre_inertial_force() { if(!current_inertial_force.is_empty()) pre_inertial_force = current_inertial_force; }
1103
1104template<sp_d T> void Factory<T>::commit_pre_displacement() { if(!current_displacement.is_empty()) pre_displacement = current_displacement; }
1105
1106template<sp_d T> void Factory<T>::commit_pre_velocity() { if(!current_velocity.is_empty()) pre_velocity = current_velocity; }
1107
1108template<sp_d T> void Factory<T>::commit_pre_acceleration() { if(!current_acceleration.is_empty()) pre_acceleration = current_acceleration; }
1109
1110template<sp_d T> void Factory<T>::commit_pre_temperature() { if(!current_temperature.is_empty()) pre_temperature = current_temperature; }
1111
1112template<sp_d T> void Factory<T>::clear_status() {
1113 access::rw(initialized) = false;
1114
1115 ninja.zeros();
1116
1117 clear_energy();
1118
1119 clear_time();
1120 clear_load_factor();
1121 clear_load();
1122 clear_settlement();
1123 clear_resistance();
1124 clear_damping_force();
1125 clear_inertial_force();
1126 clear_displacement();
1127 clear_velocity();
1128 clear_acceleration();
1129 clear_temperature();
1130 clear_auxiliary_resistance();
1131}
1132
1133template<sp_d T> void Factory<T>::clear_time() { trial_time = incre_time = current_time = 0.; }
1134
1135template<sp_d T> void Factory<T>::clear_load_factor() {
1136 if(!pre_load_factor.is_empty()) pre_load_factor.zeros();
1137 if(!trial_load_factor.is_empty()) trial_load_factor.zeros();
1138 if(!incre_load_factor.is_empty()) incre_load_factor.zeros();
1139 if(!current_load_factor.is_empty()) current_load_factor.zeros();
1140}
1141
1142template<sp_d T> void Factory<T>::clear_load() {
1143 if(!pre_load.is_empty()) pre_load.zeros();
1144 if(!trial_load.is_empty()) trial_load.zeros();
1145 if(!incre_load.is_empty()) incre_load.zeros();
1146 if(!current_load.is_empty()) current_load.zeros();
1147}
1148
1149template<sp_d T> void Factory<T>::clear_settlement() {
1150 if(!pre_settlement.is_empty()) pre_settlement.zeros();
1151 if(!trial_settlement.is_empty()) trial_settlement.zeros();
1152 if(!incre_settlement.is_empty()) incre_settlement.zeros();
1153 if(!current_settlement.is_empty()) current_settlement.zeros();
1154}
1155
1156template<sp_d T> void Factory<T>::clear_resistance() {
1157 if(!pre_resistance.is_empty()) pre_resistance.zeros();
1158 if(!trial_resistance.is_empty()) trial_resistance.zeros();
1159 if(!incre_resistance.is_empty()) incre_resistance.zeros();
1160 if(!current_resistance.is_empty()) current_resistance.zeros();
1161}
1162
1163template<sp_d T> void Factory<T>::clear_damping_force() {
1164 if(!pre_damping_force.is_empty()) pre_damping_force.zeros();
1165 if(!trial_damping_force.is_empty()) trial_damping_force.zeros();
1166 if(!incre_damping_force.is_empty()) incre_damping_force.zeros();
1167 if(!current_damping_force.is_empty()) current_damping_force.zeros();
1168}
1169
1170template<sp_d T> void Factory<T>::clear_inertial_force() {
1171 if(!pre_inertial_force.is_empty()) pre_inertial_force.zeros();
1172 if(!trial_inertial_force.is_empty()) trial_inertial_force.zeros();
1173 if(!incre_inertial_force.is_empty()) incre_inertial_force.zeros();
1174 if(!current_inertial_force.is_empty()) current_inertial_force.zeros();
1175}
1176
1177template<sp_d T> void Factory<T>::clear_displacement() {
1178 if(!pre_displacement.is_empty()) pre_displacement.zeros();
1179 if(!trial_displacement.is_empty()) trial_displacement.zeros();
1180 if(!incre_displacement.is_empty()) incre_displacement.zeros();
1181 if(!current_displacement.is_empty()) current_displacement.zeros();
1182}
1183
1184template<sp_d T> void Factory<T>::clear_velocity() {
1185 if(!pre_velocity.is_empty()) pre_velocity.zeros();
1186 if(!trial_velocity.is_empty()) trial_velocity.zeros();
1187 if(!incre_velocity.is_empty()) incre_velocity.zeros();
1188 if(!current_velocity.is_empty()) current_velocity.zeros();
1189}
1190
1191template<sp_d T> void Factory<T>::clear_acceleration() {
1192 if(!pre_acceleration.is_empty()) pre_acceleration.zeros();
1193 if(!trial_acceleration.is_empty()) trial_acceleration.zeros();
1194 if(!incre_acceleration.is_empty()) incre_acceleration.zeros();
1195 if(!current_acceleration.is_empty()) current_acceleration.zeros();
1196}
1197
1198template<sp_d T> void Factory<T>::clear_temperature() {
1199 if(!pre_temperature.is_empty()) pre_temperature.zeros();
1200 if(!trial_temperature.is_empty()) trial_temperature.zeros();
1201 if(!incre_temperature.is_empty()) incre_temperature.zeros();
1202 if(!current_temperature.is_empty()) current_temperature.zeros();
1203}
1204
1206 if(!trial_constraint_resistance.is_empty()) trial_constraint_resistance.zeros();
1207 if(!current_constraint_resistance.is_empty()) current_constraint_resistance.zeros();
1208}
1209
1210template<sp_d T> void Factory<T>::reset_status() {
1211 ninja.zeros();
1212
1213 reset_time();
1214 reset_load_factor();
1215 reset_load();
1216 reset_settlement();
1217 reset_resistance();
1218 reset_damping_force();
1219 reset_inertial_force();
1220 reset_displacement();
1221 reset_velocity();
1222 reset_acceleration();
1223 reset_temperature();
1224 reset_auxiliary_resistance();
1225}
1226
1227template<sp_d T> void Factory<T>::reset_time() {
1228 trial_time = current_time;
1229 incre_time = T(0);
1230}
1231
1232template<sp_d T> void Factory<T>::reset_load_factor() {
1233 if(trial_load_factor.is_empty()) return;
1234 trial_load_factor = current_load_factor;
1235 incre_load_factor.zeros();
1236}
1237
1238template<sp_d T> void Factory<T>::reset_load() {
1239 if(trial_load.is_empty()) return;
1240 trial_load = current_load;
1241 incre_load.zeros();
1242}
1243
1244template<sp_d T> void Factory<T>::reset_settlement() {
1245 if(trial_settlement.is_empty()) return;
1246 trial_settlement = current_settlement;
1247 incre_settlement.zeros();
1248}
1249
1250template<sp_d T> void Factory<T>::reset_resistance() {
1251 if(trial_resistance.is_empty()) return;
1252 trial_resistance = current_resistance;
1253 incre_resistance.zeros();
1254}
1255
1256template<sp_d T> void Factory<T>::reset_damping_force() {
1257 if(trial_damping_force.is_empty()) return;
1258 trial_damping_force = current_damping_force;
1259 incre_damping_force.zeros();
1260}
1261
1262template<sp_d T> void Factory<T>::reset_inertial_force() {
1263 if(trial_inertial_force.is_empty()) return;
1264 trial_inertial_force = current_inertial_force;
1265 incre_inertial_force.zeros();
1266}
1267
1268template<sp_d T> void Factory<T>::reset_displacement() {
1269 if(trial_displacement.is_empty()) return;
1270 trial_displacement = current_displacement;
1271 incre_displacement.zeros();
1272}
1273
1274template<sp_d T> void Factory<T>::reset_velocity() {
1275 if(trial_velocity.is_empty()) return;
1276 trial_velocity = current_velocity;
1277 incre_velocity.zeros();
1278}
1279
1280template<sp_d T> void Factory<T>::reset_acceleration() {
1281 if(trial_acceleration.is_empty()) return;
1282 trial_acceleration = current_acceleration;
1283 incre_acceleration.zeros();
1284}
1285
1286template<sp_d T> void Factory<T>::reset_temperature() {
1287 if(trial_temperature.is_empty()) return;
1288 trial_temperature = current_temperature;
1289 incre_temperature.zeros();
1290}
1291
1293 if(trial_constraint_resistance.is_empty()) return;
1294 trial_constraint_resistance = current_constraint_resistance;
1295}
1296
1297template<sp_d T> void Factory<T>::clear_eigen() {
1298 if(!eigenvalue.is_empty()) eigenvalue.zeros();
1299 if(!eigenvector.is_empty()) eigenvector.zeros();
1300}
1301
1302template<sp_d T> void Factory<T>::clear_mass() { if(global_mass != nullptr) global_mass->zeros(); }
1303
1304template<sp_d T> void Factory<T>::clear_damping() { if(global_damping != nullptr) global_damping->zeros(); }
1305
1306template<sp_d T> void Factory<T>::clear_stiffness() { if(global_stiffness != nullptr) global_stiffness->zeros(); }
1307
1308template<sp_d T> void Factory<T>::clear_geometry() { if(global_geometry != nullptr) global_geometry->zeros(); }
1309
1310template<sp_d T> void Factory<T>::clear_auxiliary() {
1311 n_mpc = 0;
1312 auxiliary_load.reset();
1313 auxiliary_stiffness.set_size(n_size, 0);
1314 auxiliary_resistance.reset();
1315 auxiliary_encoding.reset();
1316}
1317
1318template<sp_d T> void Factory<T>::reset() {
1319 global_mass = nullptr;
1320 global_damping = nullptr;
1321 global_stiffness = nullptr;
1322 global_geometry = nullptr;
1323}
1324
1325template<sp_d T> void Factory<T>::assemble_resistance(const Mat<T>& ER, const uvec& EI) {
1326 if(ER.is_empty()) return;
1327 for(auto I = 0llu; I < EI.n_elem; ++I) trial_resistance(EI(I)) += ER(I);
1328}
1329
1330template<sp_d T> void Factory<T>::assemble_damping_force(const Mat<T>& ER, const uvec& EI) {
1331 if(ER.is_empty()) return;
1332 for(auto I = 0llu; I < EI.n_elem; ++I) trial_damping_force(EI(I)) += ER(I);
1333}
1334
1335template<sp_d T> void Factory<T>::assemble_inertial_force(const Mat<T>& ER, const uvec& EI) {
1336 if(ER.is_empty()) return;
1337 for(auto I = 0llu; I < EI.n_elem; ++I) trial_inertial_force(EI(I)) += ER(I);
1338}
1339
1340template<sp_d T> void Factory<T>::assemble_matrix_helper(shared_ptr<MetaMat<T>>& GM, const Mat<T>& EM, const uvec& EI, const std::vector<MappingDOF>& MAP) {
1341 if(EM.is_empty()) return;
1342
1343 if(StorageScheme::BANDSYMM == storage_type || StorageScheme::SYMMPACK == storage_type) for(const auto [g_row, g_col, l_row, l_col] : MAP) GM->unsafe_at(g_row, g_col) += EM(l_row, l_col);
1344 else for(auto I = 0llu; I < EI.n_elem; ++I) for(auto J = 0llu; J < EI.n_elem; ++J) GM->unsafe_at(EI(J), EI(I)) += EM(J, I);
1345}
1346
1347template<sp_d T> void Factory<T>::assemble_mass(const Mat<T>& EM, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_mass, EM, EI, MAP); }
1348
1349template<sp_d T> void Factory<T>::assemble_damping(const Mat<T>& EC, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_damping, EC, EI, MAP); }
1350
1351template<sp_d T> void Factory<T>::assemble_stiffness(const Mat<T>& EK, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_stiffness, EK, EI, MAP); }
1352
1353template<sp_d T> void Factory<T>::assemble_geometry(const Mat<T>& EG, const uvec& EI, const std::vector<MappingDOF>& MAP) { this->assemble_matrix_helper(global_geometry, EG, EI, MAP); }
1354
1355template<sp_d T> void Factory<T>::assemble_stiffness(const SpMat<T>& EK, const uvec& EI) {
1356 if(EK.is_empty()) return;
1357 for(auto I = EK.begin(); I != EK.end(); ++I) global_stiffness->at(EI(I.row()), EI(I.col())) += *I;
1358}
1359
1360template<sp_d T> void Factory<T>::print() const {
1361 suanpan_info("A Factory object with size of {}.\n", n_size);
1362}
1363
1364template<sp_d T> unique_ptr<MetaMat<T>> Factory<T>::get_basic_container() {
1365 switch(storage_type) {
1367#ifdef SUANPAN_CUDA
1368 if(SolverType::CUDA == solver) return std::make_unique<FullMatCUDA<T>>(n_size, n_size);
1369#endif
1370 return std::make_unique<FullMat<T>>(n_size, n_size);
1372 if(SolverType::SPIKE == solver) return std::make_unique<BandMatSpike<T>>(n_size, n_lobw, n_upbw);
1373 return std::make_unique<BandMat<T>>(n_size, n_lobw, n_upbw);
1375 return std::make_unique<BandSymmMat<T>>(n_size, n_lobw);
1377 return std::make_unique<SymmPackMat<T>>(n_size);
1379 if(SolverType::MUMPS == solver) return std::make_unique<SparseMatMUMPS<T>>(n_size, n_size, n_elem);
1380 if(SolverType::LIS == solver) return std::make_unique<SparseMatLis<T>>(n_size, n_size, n_elem);
1381 if(SolverType::SUPERLU == solver) return std::make_unique<SparseMatSuperLU<T>>(n_size, n_size, n_elem);
1382#ifdef SUANPAN_MKL
1383 if(SolverType::PARDISO == solver) return std::make_unique<SparseMatPARDISO<T>>(n_size, n_size, n_elem);
1384 if(SolverType::FGMRES == solver) return std::make_unique<SparseMatFGMRES<T>>(n_size, n_size, n_elem);
1385#endif
1386#ifdef SUANPAN_CUDA
1387 if(SolverType::CUDA == solver) return std::make_unique<SparseMatCUDA<T>>(n_size, n_size, n_elem);
1388#ifdef SUANPAN_MAGMA
1389 if(SolverType::MAGMA == solver) return std::make_unique<SparseMatMAGMA<T>>(n_size, n_size, magma_setting);
1390#endif
1391#endif
1392 return std::make_unique<SparseMatSuperLU<T>>(n_size, n_size, n_elem);
1394#ifdef SUANPAN_MKL
1395 if(SolverType::FGMRES == solver) return std::make_unique<SparseSymmMatFGMRES<T>>(n_size, n_size, n_elem);
1396#endif
1397 return std::make_unique<SparseSymmMatMUMPS<T>>(n_size, n_size, n_elem);
1398 default:
1399 throw invalid_argument("need a proper storage scheme");
1400 }
1401}
1402
1403template<sp_d T> unique_ptr<MetaMat<T>> Factory<T>::get_matrix_container() {
1404 auto global_mat = get_basic_container();
1405
1406 global_mat->set_solver_setting(setting);
1407
1408 return global_mat;
1409}
1410
1411template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_mass() { return global_mass; }
1412
1413template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_damping() { return global_damping; }
1414
1415template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_stiffness() { return global_stiffness; }
1416
1417template<sp_d T> shared_ptr<MetaMat<T>>& Factory<T>::modify_geometry() { return global_geometry; }
1418
1419template<sp_d T> Col<T>& Factory<T>::modify_ninja() { return ninja; }
1420
1421template<sp_d T> Col<T>& Factory<T>::modify_sushi() { return sushi; }
1422
1423template<sp_d T> suanpan::set<uword>& Factory<T>::modify_reference_dof() { return reference_dof; }
1424
1425template<sp_d T> SpMat<T>& Factory<T>::modify_reference_load() { return reference_load; }
1426
1427template<sp_d T> uvec& Factory<T>::modify_auxiliary_encoding() { return auxiliary_encoding; }
1428
1429template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_lambda() { return auxiliary_lambda; }
1430
1431template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_resistance() { return auxiliary_resistance; }
1432
1433template<sp_d T> Col<T>& Factory<T>::modify_auxiliary_load() { return auxiliary_load; }
1434
1435template<sp_d T> SpMat<T>& Factory<T>::modify_auxiliary_stiffness() { return auxiliary_stiffness; }
1436
1437template<sp_d T> SpCol<T>& Factory<T>::modify_trial_constraint_resistance() { return trial_constraint_resistance; }
1438
1439template<sp_d T> SpCol<T>& Factory<T>::modify_current_constraint_resistance() { return current_constraint_resistance; }
1440
1441template<sp_d T> Col<T>& Factory<T>::modify_eigenvalue() { return eigenvalue; }
1442
1443template<sp_d T> Mat<T>& Factory<T>::modify_eigenvector() { return eigenvector; }
1444
1445template<sp_d T> void Factory<T>::set_trial_time(const T M) { trial_time = M; }
1446
1447template<sp_d T> void Factory<T>::set_trial_load_factor(const Col<T>& L) { trial_load_factor = L; }
1448
1449template<sp_d T> void Factory<T>::set_trial_load(const Col<T>& L) { trial_load = L; }
1450
1451template<sp_d T> void Factory<T>::set_trial_settlement(const Col<T>& S) { trial_settlement = S; }
1452
1453template<sp_d T> void Factory<T>::set_trial_resistance(const Col<T>& R) { trial_resistance = R; }
1454
1455template<sp_d T> void Factory<T>::set_trial_damping_force(const Col<T>& R) { trial_damping_force = R; }
1456
1457template<sp_d T> void Factory<T>::set_trial_inertial_force(const Col<T>& R) { trial_inertial_force = R; }
1458
1459template<sp_d T> void Factory<T>::set_trial_displacement(const Col<T>& D) { trial_displacement = D; }
1460
1461template<sp_d T> void Factory<T>::set_trial_velocity(const Col<T>& V) { trial_velocity = V; }
1462
1463template<sp_d T> void Factory<T>::set_trial_acceleration(const Col<T>& A) { trial_acceleration = A; }
1464
1465template<sp_d T> void Factory<T>::set_trial_temperature(const Col<T>& M) { trial_temperature = M; }
1466
1467template<sp_d T> void Factory<T>::set_incre_time(const T M) { incre_time = M; }
1468
1469template<sp_d T> void Factory<T>::set_incre_load_factor(const Col<T>& L) { incre_load_factor = L; }
1470
1471template<sp_d T> void Factory<T>::set_incre_load(const Col<T>& L) { incre_load = L; }
1472
1473template<sp_d T> void Factory<T>::set_incre_settlement(const Col<T>& S) { incre_settlement = S; }
1474
1475template<sp_d T> void Factory<T>::set_incre_resistance(const Col<T>& R) { incre_resistance = R; }
1476
1477template<sp_d T> void Factory<T>::set_incre_damping_force(const Col<T>& R) { incre_damping_force = R; }
1478
1479template<sp_d T> void Factory<T>::set_incre_inertial_force(const Col<T>& R) { incre_inertial_force = R; }
1480
1481template<sp_d T> void Factory<T>::set_incre_displacement(const Col<T>& D) { incre_displacement = D; }
1482
1483template<sp_d T> void Factory<T>::set_incre_velocity(const Col<T>& V) { incre_velocity = V; }
1484
1485template<sp_d T> void Factory<T>::set_incre_acceleration(const Col<T>& A) { incre_acceleration = A; }
1486
1487template<sp_d T> void Factory<T>::set_incre_temperature(const Col<T>& M) { incre_temperature = M; }
1488
1489template<sp_d T> void Factory<T>::set_current_time(const T M) { current_time = M; }
1490
1491template<sp_d T> void Factory<T>::set_current_load_factor(const Col<T>& L) { current_load_factor = L; }
1492
1493template<sp_d T> void Factory<T>::set_current_load(const Col<T>& L) { current_load = L; }
1494
1495template<sp_d T> void Factory<T>::set_current_settlement(const Col<T>& S) { current_settlement = S; }
1496
1497template<sp_d T> void Factory<T>::set_current_resistance(const Col<T>& R) { current_resistance = R; }
1498
1499template<sp_d T> void Factory<T>::set_current_damping_force(const Col<T>& R) { current_damping_force = R; }
1500
1501template<sp_d T> void Factory<T>::set_current_inertial_force(const Col<T>& R) { current_inertial_force = R; }
1502
1503template<sp_d T> void Factory<T>::set_current_displacement(const Col<T>& D) { current_displacement = D; }
1504
1505template<sp_d T> void Factory<T>::set_current_velocity(const Col<T>& V) { current_velocity = V; }
1506
1507template<sp_d T> void Factory<T>::set_current_acceleration(const Col<T>& A) { current_acceleration = A; }
1508
1509template<sp_d T> void Factory<T>::set_current_temperature(const Col<T>& M) { current_temperature = M; }
1510
1511template<sp_d T> void Factory<T>::set_pre_time(const T M) { pre_time = M; }
1512
1513template<sp_d T> void Factory<T>::set_pre_load_factor(const Col<T>& L) { pre_load_factor = L; }
1514
1515template<sp_d T> void Factory<T>::set_pre_load(const Col<T>& L) { pre_load = L; }
1516
1517template<sp_d T> void Factory<T>::set_pre_settlement(const Col<T>& S) { pre_settlement = S; }
1518
1519template<sp_d T> void Factory<T>::set_pre_resistance(const Col<T>& R) { pre_resistance = R; }
1520
1521template<sp_d T> void Factory<T>::set_pre_damping_force(const Col<T>& R) { pre_damping_force = R; }
1522
1523template<sp_d T> void Factory<T>::set_pre_inertial_force(const Col<T>& R) { pre_inertial_force = R; }
1524
1525template<sp_d T> void Factory<T>::set_pre_displacement(const Col<T>& D) { pre_displacement = D; }
1526
1527template<sp_d T> void Factory<T>::set_pre_velocity(const Col<T>& V) { pre_velocity = V; }
1528
1529template<sp_d T> void Factory<T>::set_pre_acceleration(const Col<T>& A) { pre_acceleration = A; }
1530
1531template<sp_d T> void Factory<T>::set_pre_temperature(const Col<T>& M) { pre_temperature = M; }
1532
1533template<sp_d T> T Factory<T>::get_trial_time() const { return trial_time; }
1534
1535template<sp_d T> const Col<T>& Factory<T>::get_trial_load_factor() const { return trial_load_factor; }
1536
1537template<sp_d T> const Col<T>& Factory<T>::get_trial_load() const { return trial_load; }
1538
1539template<sp_d T> const Col<T>& Factory<T>::get_trial_settlement() const { return trial_settlement; }
1540
1541template<sp_d T> const Col<T>& Factory<T>::get_trial_resistance() const { return trial_resistance; }
1542
1543template<sp_d T> const Col<T>& Factory<T>::get_trial_damping_force() const { return trial_damping_force; }
1544
1545template<sp_d T> const Col<T>& Factory<T>::get_trial_inertial_force() const { return trial_inertial_force; }
1546
1547template<sp_d T> const Col<T>& Factory<T>::get_trial_displacement() const { return trial_displacement; }
1548
1549template<sp_d T> const Col<T>& Factory<T>::get_trial_velocity() const { return trial_velocity; }
1550
1551template<sp_d T> const Col<T>& Factory<T>::get_trial_acceleration() const { return trial_acceleration; }
1552
1553template<sp_d T> const Col<T>& Factory<T>::get_trial_temperature() const { return trial_temperature; }
1554
1555template<sp_d T> T Factory<T>::get_incre_time() const { return incre_time; }
1556
1557template<sp_d T> const Col<T>& Factory<T>::get_incre_load_factor() const { return incre_load_factor; }
1558
1559template<sp_d T> const Col<T>& Factory<T>::get_incre_load() const { return incre_load; }
1560
1561template<sp_d T> const Col<T>& Factory<T>::get_incre_settlement() const { return incre_settlement; }
1562
1563template<sp_d T> const Col<T>& Factory<T>::get_incre_resistance() const { return incre_resistance; }
1564
1565template<sp_d T> const Col<T>& Factory<T>::get_incre_damping_force() const { return incre_damping_force; }
1566
1567template<sp_d T> const Col<T>& Factory<T>::get_incre_inertial_force() const { return incre_inertial_force; }
1568
1569template<sp_d T> const Col<T>& Factory<T>::get_incre_displacement() const { return incre_displacement; }
1570
1571template<sp_d T> const Col<T>& Factory<T>::get_incre_velocity() const { return incre_velocity; }
1572
1573template<sp_d T> const Col<T>& Factory<T>::get_incre_acceleration() const { return incre_acceleration; }
1574
1575template<sp_d T> const Col<T>& Factory<T>::get_incre_temperature() const { return incre_temperature; }
1576
1577template<sp_d T> T Factory<T>::get_current_time() const { return current_time; }
1578
1579template<sp_d T> const Col<T>& Factory<T>::get_current_load_factor() const { return current_load_factor; }
1580
1581template<sp_d T> const Col<T>& Factory<T>::get_current_load() const { return current_load; }
1582
1583template<sp_d T> const Col<T>& Factory<T>::get_current_settlement() const { return current_settlement; }
1584
1585template<sp_d T> const Col<T>& Factory<T>::get_current_resistance() const { return current_resistance; }
1586
1587template<sp_d T> const Col<T>& Factory<T>::get_current_damping_force() const { return current_damping_force; }
1588
1589template<sp_d T> const Col<T>& Factory<T>::get_current_inertial_force() const { return current_inertial_force; }
1590
1591template<sp_d T> const Col<T>& Factory<T>::get_current_displacement() const { return current_displacement; }
1592
1593template<sp_d T> const Col<T>& Factory<T>::get_current_velocity() const { return current_velocity; }
1594
1595template<sp_d T> const Col<T>& Factory<T>::get_current_acceleration() const { return current_acceleration; }
1596
1597template<sp_d T> const Col<T>& Factory<T>::get_current_temperature() const { return current_temperature; }
1598
1599template<sp_d T> T Factory<T>::get_pre_time() const { return pre_time; }
1600
1601template<sp_d T> const Col<T>& Factory<T>::get_pre_load_factor() const { return pre_load_factor; }
1602
1603template<sp_d T> const Col<T>& Factory<T>::get_pre_load() const { return pre_load; }
1604
1605template<sp_d T> const Col<T>& Factory<T>::get_pre_settlement() const { return pre_settlement; }
1606
1607template<sp_d T> const Col<T>& Factory<T>::get_pre_resistance() const { return pre_resistance; }
1608
1609template<sp_d T> const Col<T>& Factory<T>::get_pre_damping_force() const { return pre_damping_force; }
1610
1611template<sp_d T> const Col<T>& Factory<T>::get_pre_inertial_force() const { return pre_inertial_force; }
1612
1613template<sp_d T> const Col<T>& Factory<T>::get_pre_displacement() const { return pre_displacement; }
1614
1615template<sp_d T> const Col<T>& Factory<T>::get_pre_velocity() const { return pre_velocity; }
1616
1617template<sp_d T> const Col<T>& Factory<T>::get_pre_acceleration() const { return pre_acceleration; }
1618
1619template<sp_d T> const Col<T>& Factory<T>::get_pre_temperature() const { return pre_temperature; }
1620
1621template<sp_d T> T& Factory<T>::modify_trial_time() { return trial_time; }
1622
1623template<sp_d T> Col<T>& Factory<T>::modify_trial_load_factor() { return trial_load_factor; }
1624
1625template<sp_d T> Col<T>& Factory<T>::modify_trial_load() { return trial_load; }
1626
1627template<sp_d T> Col<T>& Factory<T>::modify_trial_settlement() { return trial_settlement; }
1628
1629template<sp_d T> Col<T>& Factory<T>::modify_trial_resistance() { return trial_resistance; }
1630
1631template<sp_d T> Col<T>& Factory<T>::modify_trial_damping_force() { return trial_damping_force; }
1632
1633template<sp_d T> Col<T>& Factory<T>::modify_trial_inertial_force() { return trial_inertial_force; }
1634
1635template<sp_d T> Col<T>& Factory<T>::modify_trial_displacement() { return trial_displacement; }
1636
1637template<sp_d T> Col<T>& Factory<T>::modify_trial_velocity() { return trial_velocity; }
1638
1639template<sp_d T> Col<T>& Factory<T>::modify_trial_acceleration() { return trial_acceleration; }
1640
1641template<sp_d T> Col<T>& Factory<T>::modify_trial_temperature() { return trial_temperature; }
1642
1643template<sp_d T> T& Factory<T>::modify_incre_time() { return incre_time; }
1644
1645template<sp_d T> Col<T>& Factory<T>::modify_incre_load_factor() { return incre_load_factor; }
1646
1647template<sp_d T> Col<T>& Factory<T>::modify_incre_load() { return incre_load; }
1648
1649template<sp_d T> Col<T>& Factory<T>::modify_incre_settlement() { return incre_settlement; }
1650
1651template<sp_d T> Col<T>& Factory<T>::modify_incre_resistance() { return incre_resistance; }
1652
1653template<sp_d T> Col<T>& Factory<T>::modify_incre_damping_force() { return incre_damping_force; }
1654
1655template<sp_d T> Col<T>& Factory<T>::modify_incre_inertial_force() { return incre_inertial_force; }
1656
1657template<sp_d T> Col<T>& Factory<T>::modify_incre_displacement() { return incre_displacement; }
1658
1659template<sp_d T> Col<T>& Factory<T>::modify_incre_velocity() { return incre_velocity; }
1660
1661template<sp_d T> Col<T>& Factory<T>::modify_incre_acceleration() { return incre_acceleration; }
1662
1663template<sp_d T> Col<T>& Factory<T>::modify_incre_temperature() { return incre_temperature; }
1664
1665template<sp_d T> T& Factory<T>::modify_current_time() { return current_time; }
1666
1667template<sp_d T> Col<T>& Factory<T>::modify_current_load_factor() { return current_load_factor; }
1668
1669template<sp_d T> Col<T>& Factory<T>::modify_current_load() { return current_load; }
1670
1671template<sp_d T> Col<T>& Factory<T>::modify_current_settlement() { return current_settlement; }
1672
1673template<sp_d T> Col<T>& Factory<T>::modify_current_resistance() { return current_resistance; }
1674
1675template<sp_d T> Col<T>& Factory<T>::modify_current_damping_force() { return current_damping_force; }
1676
1677template<sp_d T> Col<T>& Factory<T>::modify_current_inertial_force() { return current_inertial_force; }
1678
1679template<sp_d T> Col<T>& Factory<T>::modify_current_displacement() { return current_displacement; }
1680
1681template<sp_d T> Col<T>& Factory<T>::modify_current_velocity() { return current_velocity; }
1682
1683template<sp_d T> Col<T>& Factory<T>::modify_current_acceleration() { return current_acceleration; }
1684
1685template<sp_d T> Col<T>& Factory<T>::modify_current_temperature() { return current_temperature; }
1686
1687template<sp_d T> T& Factory<T>::modify_pre_time() { return pre_time; }
1688
1689template<sp_d T> Col<T>& Factory<T>::modify_pre_load_factor() { return pre_load_factor; }
1690
1691template<sp_d T> Col<T>& Factory<T>::modify_pre_load() { return pre_load; }
1692
1693template<sp_d T> Col<T>& Factory<T>::modify_pre_settlement() { return pre_settlement; }
1694
1695template<sp_d T> Col<T>& Factory<T>::modify_pre_resistance() { return pre_resistance; }
1696
1697template<sp_d T> Col<T>& Factory<T>::modify_pre_damping_force() { return pre_damping_force; }
1698
1699template<sp_d T> Col<T>& Factory<T>::modify_pre_inertial_force() { return pre_inertial_force; }
1700
1701template<sp_d T> Col<T>& Factory<T>::modify_pre_displacement() { return pre_displacement; }
1702
1703template<sp_d T> Col<T>& Factory<T>::modify_pre_velocity() { return pre_velocity; }
1704
1705template<sp_d T> Col<T>& Factory<T>::modify_pre_acceleration() { return pre_acceleration; }
1706
1707template<sp_d T> Col<T>& Factory<T>::modify_pre_temperature() { return pre_temperature; }
1708
1709template<sp_d T> void Factory<T>::update_trial_time(const T M) {
1710 trial_time = M;
1711 incre_time = trial_time - current_time;
1712}
1713
1714template<sp_d T> void Factory<T>::update_trial_load_factor(const Col<T>& L) {
1715 trial_load_factor = L;
1716 incre_load_factor = trial_load_factor - current_load_factor;
1717}
1718
1719template<sp_d T> void Factory<T>::update_trial_load(const Col<T>& L) {
1720 trial_load = L;
1721 incre_load = trial_load - current_load;
1722}
1723
1724template<sp_d T> void Factory<T>::update_trial_settlement(const Col<T>& S) {
1725 trial_settlement = S;
1726 incre_settlement = trial_settlement - current_settlement;
1727}
1728
1729template<sp_d T> void Factory<T>::update_trial_resistance(const Col<T>& R) {
1730 trial_resistance = R;
1731 incre_resistance = trial_resistance - current_resistance;
1732}
1733
1734template<sp_d T> void Factory<T>::update_trial_damping_force(const Col<T>& R) {
1735 trial_damping_force = R;
1736 incre_damping_force = trial_damping_force - current_damping_force;
1737}
1738
1739template<sp_d T> void Factory<T>::update_trial_inertial_force(const Col<T>& R) {
1740 trial_inertial_force = R;
1741 incre_inertial_force = trial_inertial_force - current_inertial_force;
1742}
1743
1744template<sp_d T> void Factory<T>::update_trial_displacement(const Col<T>& D) {
1745 trial_displacement = D;
1746 incre_displacement = trial_displacement - current_displacement;
1747}
1748
1749template<sp_d T> void Factory<T>::update_trial_velocity(const Col<T>& V) {
1750 trial_velocity = V;
1751 incre_velocity = trial_velocity - current_velocity;
1752}
1753
1754template<sp_d T> void Factory<T>::update_trial_acceleration(const Col<T>& A) {
1755 trial_acceleration = A;
1756 incre_acceleration = trial_acceleration - current_acceleration;
1757}
1758
1759template<sp_d T> void Factory<T>::update_trial_temperature(const Col<T>& M) {
1760 trial_temperature = M;
1761 incre_temperature = trial_temperature - current_temperature;
1762}
1763
1764template<sp_d T> void Factory<T>::update_incre_time(const T M) {
1765 incre_time = M;
1766 trial_time = current_time + incre_time;
1767}
1768
1769template<sp_d T> void Factory<T>::update_incre_load_factor(const Col<T>& L) {
1770 incre_load_factor = L;
1771 trial_load_factor = current_load_factor + incre_load_factor;
1772}
1773
1774template<sp_d T> void Factory<T>::update_incre_load(const Col<T>& L) {
1775 incre_load = L;
1776 trial_load = current_load + incre_load;
1777}
1778
1779template<sp_d T> void Factory<T>::update_incre_settlement(const Col<T>& S) {
1780 incre_settlement = S;
1781 trial_settlement = current_settlement + incre_settlement;
1782}
1783
1784template<sp_d T> void Factory<T>::update_incre_resistance(const Col<T>& R) {
1785 incre_resistance = R;
1786 trial_resistance = current_resistance + incre_resistance;
1787}
1788
1789template<sp_d T> void Factory<T>::update_incre_damping_force(const Col<T>& R) {
1790 incre_damping_force = R;
1791 trial_damping_force = current_damping_force + incre_damping_force;
1792}
1793
1794template<sp_d T> void Factory<T>::update_incre_inertial_force(const Col<T>& R) {
1795 incre_inertial_force = R;
1796 trial_inertial_force = current_inertial_force + incre_inertial_force;
1797}
1798
1799template<sp_d T> void Factory<T>::update_incre_displacement(const Col<T>& D) {
1800 incre_displacement = D;
1801 trial_displacement = current_displacement + incre_displacement;
1802}
1803
1804template<sp_d T> void Factory<T>::update_incre_velocity(const Col<T>& V) {
1805 incre_velocity = V;
1806 trial_velocity = current_velocity + incre_velocity;
1807}
1808
1809template<sp_d T> void Factory<T>::update_incre_acceleration(const Col<T>& A) {
1810 incre_acceleration = A;
1811 trial_acceleration = current_acceleration + incre_acceleration;
1812}
1813
1814template<sp_d T> void Factory<T>::update_incre_temperature(const Col<T>& M) {
1815 incre_temperature = M;
1816 trial_temperature = current_temperature + incre_temperature;
1817}
1818
1819template<sp_d T> void Factory<T>::update_current_time(const T M) {
1820 trial_time = current_time = M;
1821 incre_time = T(0);
1822}
1823
1824template<sp_d T> void Factory<T>::update_current_load_factor(const Col<T>& L) {
1825 trial_load_factor = current_load_factor = L;
1826 incre_load_factor.zeros();
1827}
1828
1829template<sp_d T> void Factory<T>::update_current_load(const Col<T>& L) {
1830 trial_load = current_load = L;
1831 incre_load.zeros();
1832}
1833
1834template<sp_d T> void Factory<T>::update_current_settlement(const Col<T>& S) {
1835 trial_settlement = current_settlement = S;
1836 incre_settlement.zeros();
1837}
1838
1839template<sp_d T> void Factory<T>::update_current_resistance(const Col<T>& R) {
1840 trial_resistance = current_resistance = R;
1841 incre_resistance.zeros();
1842}
1843
1844template<sp_d T> void Factory<T>::update_current_damping_force(const Col<T>& R) {
1845 trial_damping_force = current_damping_force = R;
1846 incre_damping_force.zeros();
1847}
1848
1849template<sp_d T> void Factory<T>::update_current_inertial_force(const Col<T>& R) {
1850 trial_inertial_force = current_inertial_force = R;
1851 incre_inertial_force.zeros();
1852}
1853
1854template<sp_d T> void Factory<T>::update_current_displacement(const Col<T>& D) {
1855 trial_displacement = current_displacement = D;
1856 incre_displacement.zeros();
1857}
1858
1859template<sp_d T> void Factory<T>::update_current_velocity(const Col<T>& V) {
1860 trial_velocity = current_velocity = V;
1861 incre_velocity.zeros();
1862}
1863
1864template<sp_d T> void Factory<T>::update_current_acceleration(const Col<T>& A) {
1865 trial_acceleration = current_acceleration = A;
1866 incre_acceleration.zeros();
1867}
1868
1869template<sp_d T> void Factory<T>::update_current_temperature(const Col<T>& M) {
1870 trial_temperature = current_temperature = M;
1871 incre_temperature.zeros();
1872}
1873
1874template<sp_d T> void Factory<T>::update_trial_time_by(const T M) {
1875 trial_time += M;
1876 incre_time = trial_time - current_time;
1877}
1878
1879template<sp_d T> void Factory<T>::update_trial_load_factor_by(const Col<T>& L) {
1880 trial_load_factor += L;
1881 incre_load_factor = trial_load_factor - current_load_factor;
1882}
1883
1884template<sp_d T> void Factory<T>::update_trial_load_by(const Col<T>& L) {
1885 trial_load += L;
1886 incre_load = trial_load - current_load;
1887}
1888
1889template<sp_d T> void Factory<T>::update_trial_settlement_by(const Col<T>& S) {
1890 trial_settlement += S;
1891 incre_settlement = trial_settlement - current_settlement;
1892}
1893
1894template<sp_d T> void Factory<T>::update_trial_resistance_by(const Col<T>& R) {
1895 trial_resistance += R;
1896 incre_resistance = trial_resistance - current_resistance;
1897}
1898
1899template<sp_d T> void Factory<T>::update_trial_damping_force_by(const Col<T>& R) {
1900 trial_damping_force += R;
1901 incre_damping_force = trial_damping_force - current_damping_force;
1902}
1903
1904template<sp_d T> void Factory<T>::update_trial_inertial_force_by(const Col<T>& R) {
1905 trial_inertial_force += R;
1906 incre_inertial_force = trial_inertial_force - current_inertial_force;
1907}
1908
1909template<sp_d T> void Factory<T>::update_trial_displacement_by(const Col<T>& D) {
1910 trial_displacement += D;
1911 incre_displacement = trial_displacement - current_displacement;
1912}
1913
1914template<sp_d T> void Factory<T>::update_trial_velocity_by(const Col<T>& V) {
1915 trial_velocity += V;
1916 incre_velocity = trial_velocity - current_velocity;
1917}
1918
1919template<sp_d T> void Factory<T>::update_trial_acceleration_by(const Col<T>& A) {
1920 trial_acceleration += A;
1921 incre_acceleration = trial_acceleration - current_acceleration;
1922}
1923
1924template<sp_d T> void Factory<T>::update_trial_temperature_by(const Col<T>& M) {
1925 trial_temperature += M;
1926 incre_temperature = trial_temperature - current_temperature;
1927}
1928
1929template<sp_d T> void Factory<T>::update_incre_time_by(const T M) {
1930 incre_time += M;
1931 trial_time = current_time + incre_time;
1932}
1933
1934template<sp_d T> void Factory<T>::update_incre_load_factor_by(const Col<T>& L) {
1935 incre_load_factor += L;
1936 trial_load_factor = current_load_factor + incre_load_factor;
1937}
1938
1939template<sp_d T> void Factory<T>::update_incre_load_by(const Col<T>& L) {
1940 incre_load += L;
1941 trial_load = current_load + incre_load;
1942}
1943
1944template<sp_d T> void Factory<T>::update_incre_settlement_by(const Col<T>& S) {
1945 incre_settlement += S;
1946 trial_settlement = current_settlement + incre_settlement;
1947}
1948
1949template<sp_d T> void Factory<T>::update_incre_resistance_by(const Col<T>& R) {
1950 incre_resistance += R;
1951 trial_resistance = current_resistance + incre_resistance;
1952}
1953
1954template<sp_d T> void Factory<T>::update_incre_damping_force_by(const Col<T>& R) {
1955 incre_damping_force += R;
1956 trial_damping_force = current_damping_force + incre_damping_force;
1957}
1958
1959template<sp_d T> void Factory<T>::update_incre_inertial_force_by(const Col<T>& R) {
1960 incre_inertial_force += R;
1961 trial_inertial_force = current_inertial_force + incre_inertial_force;
1962}
1963
1964template<sp_d T> void Factory<T>::update_incre_displacement_by(const Col<T>& D) {
1965 incre_displacement += D;
1966 trial_displacement = current_displacement + incre_displacement;
1967}
1968
1969template<sp_d T> void Factory<T>::update_incre_velocity_by(const Col<T>& V) {
1970 incre_velocity += V;
1971 trial_velocity = current_velocity + incre_velocity;
1972}
1973
1974template<sp_d T> void Factory<T>::update_incre_acceleration_by(const Col<T>& A) {
1975 incre_acceleration += A;
1976 trial_acceleration = current_acceleration + incre_acceleration;
1977}
1978
1979template<sp_d T> void Factory<T>::update_incre_temperature_by(const Col<T>& M) {
1980 incre_temperature += M;
1981 trial_temperature = current_temperature + incre_temperature;
1982}
1983
1984template<sp_d T> void Factory<T>::update_current_time_by(const T M) {
1985 trial_time = current_time += M;
1986 incre_time = 0.;
1987}
1988
1989template<sp_d T> void Factory<T>::update_current_load_factor_by(const Col<T>& L) {
1990 trial_load_factor = current_load_factor += L;
1991 incre_load_factor.zeros();
1992}
1993
1994template<sp_d T> void Factory<T>::update_current_load_by(const Col<T>& L) {
1995 trial_load = current_load += L;
1996 incre_load.zeros();
1997}
1998
1999template<sp_d T> void Factory<T>::update_current_settlement_by(const Col<T>& S) {
2000 trial_settlement = current_settlement += S;
2001 incre_settlement.zeros();
2002}
2003
2004template<sp_d T> void Factory<T>::update_current_resistance_by(const Col<T>& R) {
2005 trial_resistance = current_resistance += R;
2006 incre_resistance.zeros();
2007}
2008
2009template<sp_d T> void Factory<T>::update_current_damping_force_by(const Col<T>& R) {
2010 trial_damping_force = current_damping_force += R;
2011 incre_damping_force.zeros();
2012}
2013
2014template<sp_d T> void Factory<T>::update_current_inertial_force_by(const Col<T>& R) {
2015 trial_inertial_force = current_inertial_force += R;
2016 incre_inertial_force.zeros();
2017}
2018
2019template<sp_d T> void Factory<T>::update_current_displacement_by(const Col<T>& D) {
2020 trial_displacement = current_displacement += D;
2021 incre_displacement.zeros();
2022}
2023
2024template<sp_d T> void Factory<T>::update_current_velocity_by(const Col<T>& V) {
2025 trial_velocity = current_velocity += V;
2026 incre_velocity.zeros();
2027}
2028
2029template<sp_d T> void Factory<T>::update_current_acceleration_by(const Col<T>& A) {
2030 trial_acceleration = current_acceleration += A;
2031 incre_acceleration.zeros();
2032}
2033
2034template<sp_d T> void Factory<T>::update_current_temperature_by(const Col<T>& M) {
2035 trial_temperature = current_temperature += M;
2036 incre_temperature.zeros();
2037}
2038
2039#endif // FACTORY_HPP
2040
void reset(ExternalMaterialData *data, int *info)
Definition ElasticExternal.cpp:74
A Factory class.
Definition Factory.hpp:73
const bool initialized
Definition Factory.hpp:183
A MetaMat class that holds matrices.
Definition MetaMat.hpp:39
void set_pre_inertial_force(const Col< T > &)
Definition Factory.hpp:1523
T & modify_trial_time()
Definition Factory.hpp:1621
void set_eigenvalue(const Col< T > &)
Definition Factory.hpp:895
void clear_acceleration()
Definition Factory.hpp:1191
void initialize_load()
Definition Factory.hpp:795
void set_mpc(unsigned)
Definition Factory.hpp:877
shared_ptr< MetaMat< T > > & modify_stiffness()
Definition Factory.hpp:1415
Col< T > & modify_trial_inertial_force()
Definition Factory.hpp:1633
void commit_temperature()
Definition Factory.hpp:1065
bool is_sparse() const
Definition Factory.hpp:700
void set_incre_velocity(const Col< T > &)
Definition Factory.hpp:1483
uword get_entry() const
Definition Factory.hpp:666
void set_incre_temperature(const Col< T > &)
Definition Factory.hpp:1487
void assemble_stiffness(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition Factory.hpp:1351
void print() const
Definition Factory.hpp:1360
void clear_resistance()
Definition Factory.hpp:1156
void initialize_auxiliary_resistance()
Definition Factory.hpp:849
void clear_load_factor()
Definition Factory.hpp:1135
void set_incre_acceleration(const Col< T > &)
Definition Factory.hpp:1485
int initialize()
Definition Factory.hpp:734
Col< T > & modify_incre_load()
Definition Factory.hpp:1647
const Col< T > & get_current_acceleration() const
Definition Factory.hpp:1595
void update_incre_damping_force_by(const Col< T > &)
Definition Factory.hpp:1954
const Col< T > & get_pre_damping_force() const
Definition Factory.hpp:1609
void reset_displacement()
Definition Factory.hpp:1268
void update_incre_temperature_by(const Col< T > &)
Definition Factory.hpp:1979
void set_ninja(const Col< T > &)
Definition Factory.hpp:871
Col< T > & modify_trial_damping_force()
Definition Factory.hpp:1631
void update_trial_temperature_by(const Col< T > &)
Definition Factory.hpp:1924
bool is_nlgeom() const
Definition Factory.hpp:674
void update_incre_load_factor_by(const Col< T > &)
Definition Factory.hpp:1934
const Mat< T > & get_eigenvector() const
Definition Factory.hpp:963
const Col< T > & get_trial_velocity() const
Definition Factory.hpp:1549
const Col< T > & get_trial_resistance() const
Definition Factory.hpp:1541
void update_trial_load_factor_by(const Col< T > &)
Definition Factory.hpp:1879
const Col< T > & get_ninja() const
Definition Factory.hpp:899
Col< T > & modify_trial_temperature()
Definition Factory.hpp:1641
void set_pre_velocity(const Col< T > &)
Definition Factory.hpp:1527
const Col< T > & get_current_load() const
Definition Factory.hpp:1581
void set_trial_settlement(const Col< T > &)
Definition Factory.hpp:1451
void set_incre_load_factor(const Col< T > &)
Definition Factory.hpp:1469
void set_trial_time(T)
Definition Factory.hpp:1445
Col< T > & modify_incre_resistance()
Definition Factory.hpp:1651
void update_incre_load_by(const Col< T > &)
Definition Factory.hpp:1939
void set_solver_type(SolverType)
Definition Factory.hpp:676
const Col< T > & get_pre_load() const
Definition Factory.hpp:1603
void commit_acceleration()
Definition Factory.hpp:1059
AnalysisType get_analysis_type() const
Definition Factory.hpp:690
void set_current_resistance(const Col< T > &)
Definition Factory.hpp:1497
void update_incre_displacement(const Col< T > &)
Definition Factory.hpp:1799
void set_incre_displacement(const Col< T > &)
Definition Factory.hpp:1481
void update_current_acceleration_by(const Col< T > &)
Definition Factory.hpp:2029
const uvec & get_auxiliary_encoding() const
Definition Factory.hpp:907
void reset_auxiliary_resistance()
Definition Factory.hpp:1292
Col< T > & modify_pre_displacement()
Definition Factory.hpp:1701
void commit_resistance()
Definition Factory.hpp:1029
const Col< T > & get_incre_inertial_force() const
Definition Factory.hpp:1567
void set_storage_scheme(StorageScheme)
Definition Factory.hpp:692
void initialize_damping()
Definition Factory.hpp:856
void update_incre_velocity(const Col< T > &)
Definition Factory.hpp:1804
Col< T > & modify_incre_acceleration()
Definition Factory.hpp:1661
void initialize_geometry()
Definition Factory.hpp:860
void update_incre_displacement_by(const Col< T > &)
Definition Factory.hpp:1964
void update_current_damping_force_by(const Col< T > &)
Definition Factory.hpp:2009
void assemble_mass(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition Factory.hpp:1347
void set_trial_velocity(const Col< T > &)
Definition Factory.hpp:1461
void update_current_resistance(const Col< T > &)
Definition Factory.hpp:1839
void clear_inertial_force()
Definition Factory.hpp:1170
Col< T > & modify_ninja()
Definition Factory.hpp:1419
void update_trial_damping_force_by(const Col< T > &)
Definition Factory.hpp:1899
unsigned get_mpc() const
Definition Factory.hpp:903
Col< T > & modify_eigenvalue()
Definition Factory.hpp:1441
const Col< T > & get_current_settlement() const
Definition Factory.hpp:1583
const SolverSetting< double > & get_solver_setting() const
Definition Factory.hpp:682
void commit_load_factor()
Definition Factory.hpp:1011
std::mutex & get_damping_mutex()
Definition Factory.hpp:955
void update_trial_velocity_by(const Col< T > &)
Definition Factory.hpp:1914
void update_current_load_by(const Col< T > &)
Definition Factory.hpp:1994
SpCol< T > & modify_trial_constraint_resistance()
Definition Factory.hpp:1437
const shared_ptr< MetaMat< T > > & get_damping() const
Definition Factory.hpp:933
void set_pre_settlement(const Col< T > &)
Definition Factory.hpp:1517
void set_nlgeom(bool)
Definition Factory.hpp:668
StorageScheme get_storage_scheme() const
Definition Factory.hpp:698
void commit_inertial_force()
Definition Factory.hpp:1041
void update_current_time(T)
Definition Factory.hpp:1819
void initialize_resistance()
Definition Factory.hpp:807
const Col< T > & get_current_damping_force() const
Definition Factory.hpp:1587
void set_solver_setting(const SolverSetting< double > &)
Definition Factory.hpp:680
void update_trial_displacement_by(const Col< T > &)
Definition Factory.hpp:1909
void set_mass(const shared_ptr< MetaMat< T > > &)
Definition Factory.hpp:887
void reset_velocity()
Definition Factory.hpp:1274
const Col< T > & get_trial_load_factor() const
Definition Factory.hpp:1535
void set_current_load(const Col< T > &)
Definition Factory.hpp:1493
std::mutex & get_trial_load_mutex()
Definition Factory.hpp:949
void initialize_stiffness()
Definition Factory.hpp:858
void initialize_damping_force()
Definition Factory.hpp:813
T get_current_time() const
Definition Factory.hpp:1577
const Col< T > & get_pre_settlement() const
Definition Factory.hpp:1605
void update_trial_temperature(const Col< T > &)
Definition Factory.hpp:1759
void commit_velocity()
Definition Factory.hpp:1053
void set_pre_time(T)
Definition Factory.hpp:1511
T get_complementary_energy()
Definition Factory.hpp:927
void reset()
Definition Factory.hpp:1318
const Col< T > & get_incre_displacement() const
Definition Factory.hpp:1569
void set_current_damping_force(const Col< T > &)
Definition Factory.hpp:1499
void update_trial_acceleration(const Col< T > &)
Definition Factory.hpp:1754
void set_current_time(T)
Definition Factory.hpp:1489
const Col< T > & get_auxiliary_load() const
Definition Factory.hpp:913
void commit_pre_settlement()
Definition Factory.hpp:1096
Col< T > & modify_pre_settlement()
Definition Factory.hpp:1693
const shared_ptr< MetaMat< T > > & get_mass() const
Definition Factory.hpp:931
void commit_status()
Definition Factory.hpp:987
void commit_pre_load_factor()
Definition Factory.hpp:1092
const suanpan::set< uword > & get_reference_dof() const
Definition Factory.hpp:728
void set_incre_time(T)
Definition Factory.hpp:1467
shared_ptr< MetaMat< T > > & modify_damping()
Definition Factory.hpp:1413
Col< T > & modify_pre_velocity()
Definition Factory.hpp:1703
void assemble_damping(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition Factory.hpp:1349
void update_trial_inertial_force(const Col< T > &)
Definition Factory.hpp:1739
void set_current_velocity(const Col< T > &)
Definition Factory.hpp:1505
void commit_displacement()
Definition Factory.hpp:1047
void update_current_velocity_by(const Col< T > &)
Definition Factory.hpp:2024
const Col< T > & get_incre_damping_force() const
Definition Factory.hpp:1565
void update_trial_load_factor(const Col< T > &)
Definition Factory.hpp:1714
void reset_load_factor()
Definition Factory.hpp:1232
void clear_velocity()
Definition Factory.hpp:1184
void update_current_load(const Col< T > &)
Definition Factory.hpp:1829
void update_incre_time(T)
Definition Factory.hpp:1764
void update_incre_settlement_by(const Col< T > &)
Definition Factory.hpp:1944
unsigned get_reference_size() const
Definition Factory.hpp:722
void update_current_acceleration(const Col< T > &)
Definition Factory.hpp:1864
StorageScheme
Definition Factory.hpp:52
void commit_pre_temperature()
Definition Factory.hpp:1110
void commit_pre_status()
Definition Factory.hpp:1076
void update_current_time_by(T)
Definition Factory.hpp:1984
suanpan::set< uword > & modify_reference_dof()
Definition Factory.hpp:1423
void set_trial_acceleration(const Col< T > &)
Definition Factory.hpp:1463
void update_trial_damping_force(const Col< T > &)
Definition Factory.hpp:1734
const Col< T > & get_incre_settlement() const
Definition Factory.hpp:1561
void reset_inertial_force()
Definition Factory.hpp:1262
void set_reference_dof(const suanpan::set< uword > &)
Definition Factory.hpp:726
const Col< T > & get_momentum()
Definition Factory.hpp:929
void update_current_displacement(const Col< T > &)
Definition Factory.hpp:1854
void clear_settlement()
Definition Factory.hpp:1149
void set_pre_load(const Col< T > &)
Definition Factory.hpp:1515
std::mutex & get_auxiliary_stiffness_mutex()
Definition Factory.hpp:945
T get_incre_time() const
Definition Factory.hpp:1555
void commit_pre_damping_force()
Definition Factory.hpp:1100
Col< T > & modify_current_inertial_force()
Definition Factory.hpp:1677
T get_strain_energy()
Definition Factory.hpp:921
void commit_pre_inertial_force()
Definition Factory.hpp:1102
void update_trial_settlement(const Col< T > &)
Definition Factory.hpp:1724
void set_trial_inertial_force(const Col< T > &)
Definition Factory.hpp:1457
const Col< T > & get_eigenvalue() const
Definition Factory.hpp:961
void commit_pre_displacement()
Definition Factory.hpp:1104
std::mutex & get_geometry_mutex()
Definition Factory.hpp:959
void commit_pre_acceleration()
Definition Factory.hpp:1108
const Col< T > & get_incre_velocity() const
Definition Factory.hpp:1571
Col< T > & modify_current_settlement()
Definition Factory.hpp:1671
T get_trial_time() const
Definition Factory.hpp:1533
void update_current_temperature_by(const Col< T > &)
Definition Factory.hpp:2034
void set_geometry(const shared_ptr< MetaMat< T > > &)
Definition Factory.hpp:893
void set_incre_load(const Col< T > &)
Definition Factory.hpp:1471
void update_reference_size()
Definition Factory.hpp:715
void set_pre_acceleration(const Col< T > &)
Definition Factory.hpp:1529
void set_incre_inertial_force(const Col< T > &)
Definition Factory.hpp:1479
T get_kinetic_energy()
Definition Factory.hpp:923
void set_current_settlement(const Col< T > &)
Definition Factory.hpp:1495
void update_trial_settlement_by(const Col< T > &)
Definition Factory.hpp:1889
void initialize_inertial_force()
Definition Factory.hpp:819
void set_pre_resistance(const Col< T > &)
Definition Factory.hpp:1519
void commit_load()
Definition Factory.hpp:1017
void update_trial_acceleration_by(const Col< T > &)
Definition Factory.hpp:1919
const Col< T > & get_pre_acceleration() const
Definition Factory.hpp:1617
void update_incre_acceleration_by(const Col< T > &)
Definition Factory.hpp:1974
void commit_pre_load()
Definition Factory.hpp:1094
void commit_pre_velocity()
Definition Factory.hpp:1106
void update_trial_resistance_by(const Col< T > &)
Definition Factory.hpp:1894
Col< T > & modify_current_resistance()
Definition Factory.hpp:1673
AnalysisType
Definition Factory.hpp:43
void commit_pre_time()
Definition Factory.hpp:1090
shared_ptr< MetaMat< T > > & modify_geometry()
Definition Factory.hpp:1417
void commit_auxiliary_resistance()
Definition Factory.hpp:1071
void update_current_settlement(const Col< T > &)
Definition Factory.hpp:1834
void update_trial_load(const Col< T > &)
Definition Factory.hpp:1719
void update_incre_load(const Col< T > &)
Definition Factory.hpp:1774
void update_trial_time(T)
Definition Factory.hpp:1709
const shared_ptr< MetaMat< T > > & get_stiffness() const
Definition Factory.hpp:935
void commit_settlement()
Definition Factory.hpp:1023
const Col< T > & get_auxiliary_lambda() const
Definition Factory.hpp:909
void reset_load()
Definition Factory.hpp:1238
void reset_damping_force()
Definition Factory.hpp:1256
void update_current_inertial_force_by(const Col< T > &)
Definition Factory.hpp:2014
void set_trial_temperature(const Col< T > &)
Definition Factory.hpp:1465
std::mutex & get_auxiliary_encoding_mutex()
Definition Factory.hpp:939
Col< T > & modify_trial_load()
Definition Factory.hpp:1625
SolverType
Definition Factory.hpp:61
void assemble_damping_force(const Mat< T > &, const uvec &)
Definition Factory.hpp:1330
const Col< T > & get_pre_displacement() const
Definition Factory.hpp:1613
const Col< T > & get_current_inertial_force() const
Definition Factory.hpp:1589
std::mutex & get_mass_mutex()
Definition Factory.hpp:953
const Col< T > & get_pre_load_factor() const
Definition Factory.hpp:1601
Col< T > & modify_current_acceleration()
Definition Factory.hpp:1683
const Col< T > & get_trial_displacement() const
Definition Factory.hpp:1547
void update_trial_velocity(const Col< T > &)
Definition Factory.hpp:1749
void update_current_temperature(const Col< T > &)
Definition Factory.hpp:1869
void reset_temperature()
Definition Factory.hpp:1286
std::mutex & get_auxiliary_resistance_mutex()
Definition Factory.hpp:941
void set_current_acceleration(const Col< T > &)
Definition Factory.hpp:1507
void update_current_settlement_by(const Col< T > &)
Definition Factory.hpp:1999
void initialize_settlement()
Definition Factory.hpp:801
void clear_damping()
Definition Factory.hpp:1304
const Col< T > & get_incre_load_factor() const
Definition Factory.hpp:1557
void clear_energy()
Definition Factory.hpp:979
const SpMat< T > & get_reference_load() const
Definition Factory.hpp:905
void initialize_displacement()
Definition Factory.hpp:825
Col< T > & modify_auxiliary_load()
Definition Factory.hpp:1433
void reset_time()
Definition Factory.hpp:1227
const Col< T > & get_current_displacement() const
Definition Factory.hpp:1591
Col< T > & modify_current_velocity()
Definition Factory.hpp:1681
unsigned get_size() const
Definition Factory.hpp:659
T & modify_current_time()
Definition Factory.hpp:1665
void reset_settlement()
Definition Factory.hpp:1244
void clear_damping_force()
Definition Factory.hpp:1163
void set_analysis_type(AnalysisType)
Definition Factory.hpp:684
const Col< T > & get_incre_load() const
Definition Factory.hpp:1559
void set_entry(uword)
Definition Factory.hpp:661
Col< T > & modify_pre_damping_force()
Definition Factory.hpp:1697
const Col< T > & get_incre_acceleration() const
Definition Factory.hpp:1573
T get_pre_time() const
Definition Factory.hpp:1599
void clear_geometry()
Definition Factory.hpp:1308
T get_viscous_energy()
Definition Factory.hpp:925
std::mutex & get_trial_constraint_resistance_mutex()
Definition Factory.hpp:947
const Col< T > & get_current_temperature() const
Definition Factory.hpp:1597
void update_incre_time_by(T)
Definition Factory.hpp:1929
void set_reference_size(unsigned)
Definition Factory.hpp:717
Col< T > & modify_trial_settlement()
Definition Factory.hpp:1627
void initialize_temperature()
Definition Factory.hpp:843
void set_current_displacement(const Col< T > &)
Definition Factory.hpp:1503
void update_reference_dof(const uvec &)
Definition Factory.hpp:724
const Col< T > & get_trial_inertial_force() const
Definition Factory.hpp:1545
void set_current_inertial_force(const Col< T > &)
Definition Factory.hpp:1501
SpCol< T > & modify_current_constraint_resistance()
Definition Factory.hpp:1439
Col< T > & modify_incre_load_factor()
Definition Factory.hpp:1645
void set_trial_load(const Col< T > &)
Definition Factory.hpp:1449
Col< T > & modify_auxiliary_resistance()
Definition Factory.hpp:1431
Col< T > & modify_pre_resistance()
Definition Factory.hpp:1695
Col< T > & modify_incre_temperature()
Definition Factory.hpp:1663
T & modify_pre_time()
Definition Factory.hpp:1687
void update_incre_temperature(const Col< T > &)
Definition Factory.hpp:1814
void update_current_velocity(const Col< T > &)
Definition Factory.hpp:1859
void set_error(T)
Definition Factory.hpp:730
void assemble_geometry(const Mat< T > &, const uvec &, const std::vector< MappingDOF > &)
Definition Factory.hpp:1353
const Col< T > & get_trial_settlement() const
Definition Factory.hpp:1539
void set_eigenvector(const Mat< T > &)
Definition Factory.hpp:897
void initialize_load_factor()
Definition Factory.hpp:785
Col< T > & modify_incre_inertial_force()
Definition Factory.hpp:1655
const Col< T > & get_trial_load() const
Definition Factory.hpp:1537
uvec & modify_auxiliary_encoding()
Definition Factory.hpp:1427
void update_current_displacement_by(const Col< T > &)
Definition Factory.hpp:2019
void set_trial_resistance(const Col< T > &)
Definition Factory.hpp:1453
void set_pre_load_factor(const Col< T > &)
Definition Factory.hpp:1513
void update_current_damping_force(const Col< T > &)
Definition Factory.hpp:1844
Col< T > & modify_current_damping_force()
Definition Factory.hpp:1675
void update_current_load_factor(const Col< T > &)
Definition Factory.hpp:1824
void commit_energy()
Definition Factory.hpp:965
void commit_damping_force()
Definition Factory.hpp:1035
T get_error() const
Definition Factory.hpp:732
void update_trial_time_by(T)
Definition Factory.hpp:1874
void update_incre_acceleration(const Col< T > &)
Definition Factory.hpp:1809
const SpMat< T > & get_auxiliary_stiffness() const
Definition Factory.hpp:915
void assemble_inertial_force(const Mat< T > &, const uvec &)
Definition Factory.hpp:1335
void clear_status()
Definition Factory.hpp:1112
void set_bandwidth(unsigned, unsigned)
Definition Factory.hpp:702
void update_incre_velocity_by(const Col< T > &)
Definition Factory.hpp:1969
const Col< T > & get_current_load_factor() const
Definition Factory.hpp:1579
void commit_time()
Definition Factory.hpp:1006
void clear_mass()
Definition Factory.hpp:1302
void get_bandwidth(unsigned &, unsigned &) const
Definition Factory.hpp:710
const Col< T > & get_pre_inertial_force() const
Definition Factory.hpp:1611
void commit_pre_resistance()
Definition Factory.hpp:1098
void reset_resistance()
Definition Factory.hpp:1250
SolverType get_solver_type() const
Definition Factory.hpp:678
Col< T > & modify_pre_temperature()
Definition Factory.hpp:1707
std::mutex & get_auxiliary_load_mutex()
Definition Factory.hpp:943
const Col< T > & get_trial_acceleration() const
Definition Factory.hpp:1551
void assemble_resistance(const Mat< T > &, const uvec &)
Definition Factory.hpp:1325
void set_incre_damping_force(const Col< T > &)
Definition Factory.hpp:1477
Col< T > & modify_incre_displacement()
Definition Factory.hpp:1657
Col< T > & modify_pre_load()
Definition Factory.hpp:1691
const SpCol< T > & get_trial_constraint_resistance() const
Definition Factory.hpp:917
void update_incre_damping_force(const Col< T > &)
Definition Factory.hpp:1789
const Col< T > & get_pre_velocity() const
Definition Factory.hpp:1615
const Col< T > & get_current_resistance() const
Definition Factory.hpp:1585
void initialize_velocity()
Definition Factory.hpp:831
std::mutex & get_stiffness_mutex()
Definition Factory.hpp:957
Col< T > & modify_pre_acceleration()
Definition Factory.hpp:1705
void update_incre_inertial_force_by(const Col< T > &)
Definition Factory.hpp:1959
const Col< T > & get_current_velocity() const
Definition Factory.hpp:1593
void update_incre_load_factor(const Col< T > &)
Definition Factory.hpp:1769
SpMat< T > & modify_auxiliary_stiffness()
Definition Factory.hpp:1435
void set_incre_resistance(const Col< T > &)
Definition Factory.hpp:1475
Col< T > & modify_trial_velocity()
Definition Factory.hpp:1637
Col< T > & modify_trial_displacement()
Definition Factory.hpp:1635
const Col< T > & get_pre_temperature() const
Definition Factory.hpp:1619
Col< T > & modify_trial_acceleration()
Definition Factory.hpp:1639
void clear_load()
Definition Factory.hpp:1142
const Col< T > & get_incre_temperature() const
Definition Factory.hpp:1575
Col< T > & modify_pre_load_factor()
Definition Factory.hpp:1689
Factory(unsigned=0, AnalysisType=AnalysisType::NONE, StorageScheme=StorageScheme::FULL)
Definition Factory.hpp:648
Col< T > & modify_trial_load_factor()
Definition Factory.hpp:1623
void set_size(unsigned)
Definition Factory.hpp:653
void update_trial_displacement(const Col< T > &)
Definition Factory.hpp:1744
const Col< T > & get_incre_resistance() const
Definition Factory.hpp:1563
Col< T > & modify_trial_resistance()
Definition Factory.hpp:1629
const Col< T > & get_sushi() const
Definition Factory.hpp:901
const Col< T > & get_auxiliary_resistance() const
Definition Factory.hpp:911
const Col< T > & get_trial_damping_force() const
Definition Factory.hpp:1543
Mat< T > & modify_eigenvector()
Definition Factory.hpp:1443
void update_current_resistance_by(const Col< T > &)
Definition Factory.hpp:2004
void set_pre_temperature(const Col< T > &)
Definition Factory.hpp:1531
void update_trial_load_by(const Col< T > &)
Definition Factory.hpp:1884
Col< T > & modify_auxiliary_lambda()
Definition Factory.hpp:1429
const SpCol< T > & get_current_constraint_resistance() const
Definition Factory.hpp:919
void clear_auxiliary_resistance()
Definition Factory.hpp:1205
void update_incre_resistance_by(const Col< T > &)
Definition Factory.hpp:1949
void update_sushi_by(const Col< T > &)
Definition Factory.hpp:875
void set_incre_settlement(const Col< T > &)
Definition Factory.hpp:1473
Col< T > & modify_pre_inertial_force()
Definition Factory.hpp:1699
void set_sushi(const Col< T > &)
Definition Factory.hpp:873
shared_ptr< MetaMat< T > > & modify_mass()
Definition Factory.hpp:1411
void update_incre_settlement(const Col< T > &)
Definition Factory.hpp:1779
void set_pre_displacement(const Col< T > &)
Definition Factory.hpp:1525
void set_current_load_factor(const Col< T > &)
Definition Factory.hpp:1491
void set_trial_displacement(const Col< T > &)
Definition Factory.hpp:1459
void update_trial_resistance(const Col< T > &)
Definition Factory.hpp:1729
Col< T > & modify_current_load_factor()
Definition Factory.hpp:1667
void clear_eigen()
Definition Factory.hpp:1297
void clear_auxiliary()
Definition Factory.hpp:1310
void set_reference_load(const SpMat< T > &)
Definition Factory.hpp:885
T & modify_incre_time()
Definition Factory.hpp:1643
void update_incre_resistance(const Col< T > &)
Definition Factory.hpp:1784
SpMat< T > & modify_reference_load()
Definition Factory.hpp:1425
Col< T > & modify_incre_damping_force()
Definition Factory.hpp:1653
Col< T > & modify_incre_settlement()
Definition Factory.hpp:1649
void update_current_load_factor_by(const Col< T > &)
Definition Factory.hpp:1989
void set_trial_load_factor(const Col< T > &)
Definition Factory.hpp:1447
void initialize_eigen()
Definition Factory.hpp:866
void update_incre_inertial_force(const Col< T > &)
Definition Factory.hpp:1794
void clear_time()
Definition Factory.hpp:1133
void set_stiffness(const shared_ptr< MetaMat< T > > &)
Definition Factory.hpp:891
void clear_displacement()
Definition Factory.hpp:1177
const shared_ptr< MetaMat< T > > & get_geometry() const
Definition Factory.hpp:937
Col< T > & modify_current_load()
Definition Factory.hpp:1669
void reset_status()
Definition Factory.hpp:1210
Col< T > & modify_current_temperature()
Definition Factory.hpp:1685
void update_current_inertial_force(const Col< T > &)
Definition Factory.hpp:1849
const Col< T > & get_trial_temperature() const
Definition Factory.hpp:1553
void set_trial_damping_force(const Col< T > &)
Definition Factory.hpp:1455
std::mutex & get_trial_settlement_mutex()
Definition Factory.hpp:951
void initialize_mass()
Definition Factory.hpp:854
Col< T > & modify_sushi()
Definition Factory.hpp:1421
void set_pre_damping_force(const Col< T > &)
Definition Factory.hpp:1521
void clear_temperature()
Definition Factory.hpp:1198
Col< T > & modify_current_displacement()
Definition Factory.hpp:1679
void initialize_acceleration()
Definition Factory.hpp:837
void update_trial_inertial_force_by(const Col< T > &)
Definition Factory.hpp:1904
void set_current_temperature(const Col< T > &)
Definition Factory.hpp:1509
Col< T > & modify_incre_velocity()
Definition Factory.hpp:1659
void set_damping(const shared_ptr< MetaMat< T > > &)
Definition Factory.hpp:889
void reset_acceleration()
Definition Factory.hpp:1280
const Col< T > & get_pre_resistance() const
Definition Factory.hpp:1607
void clear_stiffness()
Definition Factory.hpp:1306
std::set< T > set
Definition container.h:54
Definition SolverSetting.hpp:40
#define suanpan_info
Definition suanPan.h:305