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