KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
explicit_strategy.h
Go to the documentation of this file.
1 /*
2 ==============================================================================
3 KratosStructuralApplication
4 A library based on:
5 Kratos
6 A General Purpose Software for Multi-Physics Finite Element Analysis
7 Version 1.0 (Released on march 05, 2007).
8 
9 Copyright 2007
10 Pooyan Dadvand, Riccardo Rossi, Janosch Stascheit, Felix Nagel
11 pooyan@cimne.upc.edu
12 rrossi@cimne.upc.edu
13 janosch.stascheit@rub.de
14 nagel@sd.rub.de
15 - CIMNE (International Center for Numerical Methods in Engineering),
16 Gran Capita' s/n, 08034 Barcelona, Spain
17 - Ruhr-University Bochum, Institute for Structural Mechanics, Germany
18 
19 
20 Permission is hereby granted, free of charge, to any person obtaining
21 a copy of this software and associated documentation files (the
22 "Software"), to deal in the Software without restriction, including
23 without limitation the rights to use, copy, modify, merge, publish,
24 distribute, sublicense and/or sell copies of the Software, and to
25 permit persons to whom the Software is furnished to do so, subject to
26 the following condition:
27 
28 Distribution of this code for any commercial purpose is permissible
29 ONLY BY DIRECT ARRANGEMENT WITH THE COPYRIGHT OWNERS.
30 
31 The above copyright notice and this permission notice shall be
32 included in all copies or substantial portions of the Software.
33 
34 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
35 EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
36 MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
37 IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
38 CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
39 TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
40 SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
41 
42 ==============================================================================
43 */
44 /* *********************************************************
45 *
46 * Last Modified by: $Author: Nelson $
47 * Date: $Date: 2009-09-18 $
48 * Revision: $Revision: 1.0 $
49 *
50 * ***********************************************************/
51 
52 
53 
54 #if !defined(PFEM2_EXPLICIT_STRATEGY)
55 #define KRATOS_PFEM2_EXPLICIT_STRATEGY
56 
57 
58 /* System includes */
59 #include <string>
60 #include <iostream>
61 #include <algorithm>
62 
64 
65 /* External includes */
66 #ifdef _OPENMP
67 #include <omp.h>
68 #endif
69 
70 #include "boost/smart_ptr.hpp"
71 
72 
73 /* Project includes */
74 #include "includes/define.h"
75 #include "includes/model_part.h"
76 #include "includes/element.h"
79 //#include "solving_strategies/builder_and_solvers/residualbased_elimination_builder_and_solver.h"
80 #include "includes/variables.h"
81 #include "includes/cfd_variables.h"
82 #include "containers/array_1d.h"
84 //#include "custom_utilities/neighbours_calculator.h"
85 //#include "custom_elements/2fluid_2d.h"
86 //#include "custom_elements/2fluid_3d.h"
87 
88 namespace Kratos
89 {
90 
91  template<
92  class TSparseSpace,
93  class TDenseSpace,
94  class TLinearSolver>
95  class PFEM2_Explicit_Strategy : public ImplicitSolvingStrategy<TSparseSpace,TDenseSpace,TLinearSolver>
96  {
97 
98  public:
99 
101 
103 
104  typedef typename BaseType::TDataType TDataType;
105 
106  typedef TSparseSpace SparseSpaceType;
107 
109 
111 
113 
115 
117 
119 
121 
123 
125 
127 
129 
131 
132  typedef ConditionsContainerType::iterator ConditionsContainerIterator;
133 
135 
137 
139 
140 
141  //typedef Element::Pointer ParticlePointer;
142  //typedef typename std::vector<ParticlePointer> ParticlePointerVector;
143  //typedef typename std::vector<ParticlePointer>::iterator ParticlePointerIterator;
144 
145  //typedef GlobalPointersVector<Element > ParticleWeakVector;
146  //typedef GlobalPointersVector<Element >::iterator ParticleWeakIterator;
147 
148 
149 
150 
153  const int dimension,
154  // const int damp_type,
155  //const double damping_ratio,
156  // const bool virtual_mass,
157  // const double contact_stiffness_ratio,
158  // const double max_delta_time,
159  //const bool CalculateReactions,
160  //const bool ComputeFemFemContact,
161  const bool MoveMeshFlag
162  //typename TLinearSolver::Pointer pNewLinearSolver,
163  //typename TSchemeType::Pointer pScheme,
164  //typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver
165  )
166 
167  : ImplicitSolvingStrategy<TSparseSpace,TDenseSpace,TLinearSolver>(model_part, MoveMeshFlag)
168  {
169 
170  }
171 
173 
174 
175 
176  //********************************************
177  //********************************************
178  inline void CreatePartition(unsigned int number_of_threads, const int number_of_rows, vector<unsigned int>& partitions)
179  {
180  partitions.resize(number_of_threads+1);
181  int partition_size = number_of_rows / number_of_threads;
182  partitions[0] = 0;
183  partitions[number_of_threads] = number_of_rows;
184  for(unsigned int i = 1; i<number_of_threads; i++)
185  partitions[i] = partitions[i-1] + partition_size ;
186  }
187 
188 
190 {
191 
192  KRATOS_TRY
193  ModelPart& r_model_part = BaseType::GetModelPart();
194  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
195  ElementsArrayType& pElements = r_model_part.Elements();
196 
197 
198  typename ElementsArrayType::iterator it_begin = pElements.ptr_begin() ;
199  typename ElementsArrayType::iterator it_end = pElements.ptr_end();
200  for (ElementsArrayType::iterator it = it_begin; it != it_end; ++it)
201  {
202  it->AddExplicitContribution(CurrentProcessInfo);
203  }
204 
205 
206  KRATOS_CATCH("")
207 }
208 
209  //SPECIFIC FUNCTIONS FOR MY APPLICATION
210 void InitializeSolutionStep() override
211 {
212  KRATOS_TRY
213 
214  ModelPart& r_model_part = BaseType::GetModelPart();
215  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
216 
217  switch ( CurrentProcessInfo[FRACTIONAL_STEP] )
218  {
219  case 0:
220  {
221  SetToZeroVariablesInViscousIterations(CurrentProcessInfo);
222  break;
223  }
224  case 3:
225  {
226  SetToZeroVariablesInPresureIterations(CurrentProcessInfo);
227  break;
228  }
229  case 4:
230  {
231  SetToZeroVariablesForVolumetricStrain(CurrentProcessInfo);
232  break;
233  }
234  case 5:
235  {
236  SetToZeroVariablesForPressure(CurrentProcessInfo);
237  break;
238  }
239  case 6:
240  {
242  break;
243  }
244  case 7:
245  {
246  SetToZeroMassAndArea(CurrentProcessInfo);
247  break;
248  }
249  case 10:
250  {
251  SetToZeroVariablesInPresureProjection(CurrentProcessInfo);
252  break;
253  }
254 
255 
256 
257  default:
258  {
259  KRATOS_THROW_ERROR(std::logic_error,"Unexpected value for FRACTIONAL_STEP index: ", CurrentProcessInfo[FRACTIONAL_STEP]);
260  }
261  }
262  KRATOS_CATCH("")
263 }
264 
265 void FinalizeSolutionStep() override
266 {
267  KRATOS_TRY
268 
269  ModelPart& r_model_part = BaseType::GetModelPart();
270  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
271 
272  switch ( CurrentProcessInfo[FRACTIONAL_STEP] )
273  {
274  case 0:
275  {
277  break;
278  }
279  case 3:
280  {
282  break;
283  }
284  case 4:
285  {
286  UpdateLoopForVolumetricStrain(CurrentProcessInfo);
287  break;
288  }
289  case 5:
290  {
291  UpdateLoopForPressure(CurrentProcessInfo);
292  break;
293  }
294  case 6:
295  {
296  UpdateLoopForPressureViscousCorrection(CurrentProcessInfo);
297  break;
298  }
299  case 7:
300  {
301  UpdateLoopForMassAndArea(CurrentProcessInfo);
302  break;
303  }
304  case 10:
305  {
306  NormalizePressureProjection(CurrentProcessInfo);
307  break;
308  }
309 
310  default:
311  {
312  KRATOS_THROW_ERROR(std::logic_error,"Unexpected value for FRACTIONAL_STEP index: ", CurrentProcessInfo[FRACTIONAL_STEP]);
313  }
314  }
315  KRATOS_CATCH("")
316 }
317 
318 
319 
320 
321 //VISCOUS ITERATIONS
322 
324 {
325  KRATOS_TRY
326 
327  ModelPart& r_model_part = BaseType::GetModelPart();
328  NodesArrayType& pNodes = r_model_part.Nodes();
329 
330  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
331  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
332 
333  #ifdef _OPENMP
334  int number_of_threads = omp_get_max_threads();
335  #else
336  int number_of_threads = 1;
337  #endif
338 
339  vector<unsigned int> node_partition;
340  CreatePartition(number_of_threads, pNodes.size(), node_partition);
341 
342  #pragma omp parallel for
343  for(int k=0; k<number_of_threads; k++)
344  {
345  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
346  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
347 
348  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
349  {
350  noalias(i->FastGetSolutionStepValue(RHS)) = ZeroVector(3);
351  i->FastGetSolutionStepValue(NODAL_MASS)=0.0;
352  }
353  }
354 
355  KRATOS_CATCH("")
356 }
357 
359 {
360  KRATOS_TRY
361 
362  ModelPart& r_model_part = BaseType::GetModelPart();
363  NodesArrayType& pNodes = r_model_part.Nodes();
364  //const double delta_t = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
365  //const double factor = delta_t;
366 
367  #ifdef _OPENMP
368  int number_of_threads = omp_get_max_threads();
369  #else
370  int number_of_threads = 1;
371  #endif
372 
373  vector<unsigned int> node_partition;
374  CreatePartition(number_of_threads, pNodes.size(), node_partition);
375 
376  #pragma omp parallel for
377  for(int k=0; k<number_of_threads; k++)
378  {
379  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
380  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
381 
382  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
383  {
384  array_1d<double,3>& rhs = (i->FastGetSolutionStepValue(RHS));
385  array_1d<double,3>& node_update_variable = (i->FastGetSolutionStepValue(VELOCITY));
386  //noalias(node_update_variable) += factor* acceleration ;
387  noalias(node_update_variable) = rhs /(i->FastGetSolutionStepValue(NODAL_MASS)) ;
388 
389  }
390  }
391 
392  KRATOS_CATCH("")
393 }
394 
395 
396 //PRESSURE ITERATIONS
397 
399 {
400  KRATOS_TRY
401 
402  ModelPart& r_model_part = BaseType::GetModelPart();
403  NodesArrayType& pNodes = r_model_part.Nodes();
404 
405  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
406  const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
407 
408  #ifdef _OPENMP
409  int number_of_threads = omp_get_max_threads();
410  #else
411  int number_of_threads = 1;
412  #endif
413 
414  vector<unsigned int> node_partition;
415  CreatePartition(number_of_threads, pNodes.size(), node_partition);
416 
417  #pragma omp parallel for
418  for(int k=0; k<number_of_threads; k++)
419  {
420  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
421  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
422 
423  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
424  {
425 
426  if (iteration_number == 1)
427  noalias(i->FastGetSolutionStepValue(ACCELERATION)) = ZeroVector(3);
428  else// if (iteration_number==2)
429  noalias(i->FastGetSolutionStepValue(ACCELERATION)) = i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO); //second order cos we are in the second (or higher) iteration.
430  //else if (iteration_number==3)
431  // noalias(inode->FastGetSolutionStepValue(ACCELERATION)) = - inode->FastGetSolutionStepValue(ACCELERATION) + inode->FastGetSolutionStepValue(PRESS_PROJ_NO_RO);
432  noalias(i->FastGetSolutionStepValue(PRESS_PROJ)) = ZeroVector(3);
433 
434  noalias(i->FastGetSolutionStepValue(PRESS_PROJ)) = ZeroVector(3);
435  noalias(i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO)) = ZeroVector(3);
436 
437 
438  //noalias(in->GetSolutionStepValue(PRESSURE,1))=in->FastGetSolutionStepValue(PRESSURE);
439  //noalias(in->FastGetSolutionStepValue(PRESSURE)) = 0.0;
440  i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
441  i->FastGetSolutionStepValue(NODAL_MASS)=0.0;
442  //i->FastGetSolutionStepValue(MASS)=ZeroVector(3);
443 
444  //i->FastGetSolutionStepValue(PRESSURE)=0.0;
445  i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN)=0.0;
446 
447  }
448  }
449 
450  KRATOS_CATCH("")
451 }
452 
453 
455 {
456  KRATOS_TRY
457 
458  ModelPart& r_model_part = BaseType::GetModelPart();
459  NodesArrayType& pNodes = r_model_part.Nodes();
460 
461  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
462  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
463 
464  #ifdef _OPENMP
465  int number_of_threads = omp_get_max_threads();
466  #else
467  int number_of_threads = 1;
468  #endif
469 
470  vector<unsigned int> node_partition;
471  CreatePartition(number_of_threads, pNodes.size(), node_partition);
472 
473  #pragma omp parallel for
474  for(int k=0; k<number_of_threads; k++)
475  {
476  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
477  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
478 
479  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
480  {
481  i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
482  i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN)=0.0;
483  }
484  }
485 
486  KRATOS_CATCH("")
487 }
488 
490 {
491  KRATOS_TRY
492 
493  ModelPart& r_model_part = BaseType::GetModelPart();
494  NodesArrayType& pNodes = r_model_part.Nodes();
495 
496  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
497  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
498 
499  #ifdef _OPENMP
500  int number_of_threads = omp_get_max_threads();
501  #else
502  int number_of_threads = 1;
503  #endif
504 
505  vector<unsigned int> node_partition;
506  CreatePartition(number_of_threads, pNodes.size(), node_partition);
507 
508  #pragma omp parallel for
509  for(int k=0; k<number_of_threads; k++)
510  {
511  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
512  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
513 
514  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
515  {
516  i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
517  i->FastGetSolutionStepValue(PRESSURE)=0.0;
518  }
519  }
520 
521  KRATOS_CATCH("")
522 }
523 
525 {
526  KRATOS_TRY
527 
528  ModelPart& r_model_part = BaseType::GetModelPart();
529  NodesArrayType& pNodes = r_model_part.Nodes();
530  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
531  #ifdef _OPENMP
532  int number_of_threads = omp_get_max_threads();
533  #else
534  int number_of_threads = 1;
535  #endif
536 
537  vector<unsigned int> node_partition;
538  CreatePartition(number_of_threads, pNodes.size(), node_partition);
539 
540  #pragma omp parallel for
541  for(int k=0; k<number_of_threads; k++)
542  {
543  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
544  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
545 
546  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
547  {
548  //normalizing variables:
549  array_1d<double,3>& press_proj_no_ro = (i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO));
550  array_1d<double,3>& press_proj_stabilization = (i->FastGetSolutionStepValue(PRESS_PROJ));
551  double& mass = (i->FastGetSolutionStepValue(NODAL_MASS)); //this already includes 1/delta_t
552  double& area = (i->FastGetSolutionStepValue(NODAL_AREA));
553  press_proj_no_ro /= mass; //so this is already (pres_proj / mass) * delta_t;
554  press_proj_stabilization /= area;
555 
556 
557 
558  //updating acceleration
559  array_1d<double,3>& acceleration = (i->FastGetSolutionStepValue(ACCELERATION));
560  noalias(acceleration) -= (press_proj_no_ro);
561  //updating variable
562  array_1d<double,3>& velocity = (i->FastGetSolutionStepValue(VELOCITY));
563  noalias(velocity) += (acceleration) ; //the nodal mass includes the delta_time, so this is actually acceleration*delta_t
564 
565  i->FastGetSolutionStepValue(PREVIOUS_ITERATION_PRESSURE)=i->FastGetSolutionStepValue(PRESSURE);
566  }
567  }
568 
569  KRATOS_CATCH("")
570 }
571 
572 
573 //to calculate only the pressure projection:
574 
576 {
577  KRATOS_TRY
578 
579  ModelPart& r_model_part = BaseType::GetModelPart();
580  NodesArrayType& pNodes = r_model_part.Nodes();
581 
582  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
583  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
584 
585  #ifdef _OPENMP
586  int number_of_threads = omp_get_max_threads();
587  #else
588  int number_of_threads = 1;
589  #endif
590 
591  vector<unsigned int> node_partition;
592  CreatePartition(number_of_threads, pNodes.size(), node_partition);
593 
594  #pragma omp parallel for
595  for(int k=0; k<number_of_threads; k++)
596  {
597  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
598  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
599 
600  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
601  {
602 
603  noalias(i->FastGetSolutionStepValue(PRESS_PROJ)) = ZeroVector(3);
604 
605  i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
606 
607  }
608  }
609 
610  KRATOS_CATCH("")
611 }
612 
613 
614 void NormalizePressureProjection(ProcessInfo& CurrentProcessInfo)
615 {
616  KRATOS_TRY
617 
618  ModelPart& r_model_part = BaseType::GetModelPart();
619  NodesArrayType& pNodes = r_model_part.Nodes();
620  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
621  #ifdef _OPENMP
622  int number_of_threads = omp_get_max_threads();
623  #else
624  int number_of_threads = 1;
625  #endif
626 
627  vector<unsigned int> node_partition;
628  CreatePartition(number_of_threads, pNodes.size(), node_partition);
629 
630  #pragma omp parallel for
631  for(int k=0; k<number_of_threads; k++)
632  {
633  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
634  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
635 
636  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
637  {
638  //normalizing variables:
639  //array_1d<double,3>& press_proj_no_ro = (i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO));
640  array_1d<double,3>& press_proj_stabilization = (i->FastGetSolutionStepValue(PRESS_PROJ));
641  //double& mass = (i->FastGetSolutionStepValue(NODAL_MASS));
642  //array_1d<double,3>& vectorial_mass=(i->FastGetSolutionStepValue(MASS));
643  double& area = (i->FastGetSolutionStepValue(NODAL_AREA));
644  //press_proj_no_ro /= mass;
645  //press_proj_no_ro(0) /= vectorial_mass(0);
646  //press_proj_no_ro(1) /= vectorial_mass(1);
647  //press_proj_no_ro(2) /= vectorial_mass(2)+1e-20;
648  press_proj_stabilization /= area;
649 
650  //(i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN))=i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN)/area;
651  //(i->FastGetSolutionStepValue(PRESSURE))=i->FastGetSolutionStepValue(PRESSURE)/area;
652  //(i->FastGetSolutionStepValue(ELASTIC_PRESSURE))=i->FastGetSolutionStepValue(ELASTIC_PRESSURE)/area;
653 
654  }
655  }
656 
657  KRATOS_CATCH("")
658 }
659 
661 {
662  KRATOS_TRY
663 
664  ModelPart& r_model_part = BaseType::GetModelPart();
665  NodesArrayType& pNodes = r_model_part.Nodes();
666  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
667  #ifdef _OPENMP
668  int number_of_threads = omp_get_max_threads();
669  #else
670  int number_of_threads = 1;
671  #endif
672 
673  vector<unsigned int> node_partition;
674  CreatePartition(number_of_threads, pNodes.size(), node_partition);
675 
676  #pragma omp parallel for
677  for(int k=0; k<number_of_threads; k++)
678  {
679  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
680  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
681 
682  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
683  {
684  double& area = (i->FastGetSolutionStepValue(NODAL_AREA));
685  (i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN))=i->FastGetSolutionStepValue(VOLUMETRIC_STRAIN)/area;
686  }
687  }
688 
689  KRATOS_CATCH("")
690 }
691 
692 void UpdateLoopForPressure(ProcessInfo& CurrentProcessInfo)
693 {
694  KRATOS_TRY
695 
696  ModelPart& r_model_part = BaseType::GetModelPart();
697  NodesArrayType& pNodes = r_model_part.Nodes();
698  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
699  #ifdef _OPENMP
700  int number_of_threads = omp_get_max_threads();
701  #else
702  int number_of_threads = 1;
703  #endif
704 
705  vector<unsigned int> node_partition;
706  CreatePartition(number_of_threads, pNodes.size(), node_partition);
707 
708  #pragma omp parallel for
709  for(int k=0; k<number_of_threads; k++)
710  {
711  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
712  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
713 
714  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
715  {
716  double& area = (i->FastGetSolutionStepValue(NODAL_AREA));
717  (i->FastGetSolutionStepValue(PRESSURE))=i->FastGetSolutionStepValue(PRESSURE)/area;
718  }
719  }
720 
721  KRATOS_CATCH("")
722 }
723 
724 
725 
726 
727 //PRESSURE VISCOUS CORRECTION
728 
730 {
731  KRATOS_TRY
732 
733  ModelPart& r_model_part = BaseType::GetModelPart();
734  NodesArrayType& pNodes = r_model_part.Nodes();
735 
736  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
737  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
738 
739  #ifdef _OPENMP
740  int number_of_threads = omp_get_max_threads();
741  #else
742  int number_of_threads = 1;
743  #endif
744 
745  vector<unsigned int> node_partition;
746  CreatePartition(number_of_threads, pNodes.size(), node_partition);
747 
748  #pragma omp parallel for
749  for(int k=0; k<number_of_threads; k++)
750  {
751  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
752  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
753 
754  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
755  {
756  i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
757  i->FastGetSolutionStepValue(NODAL_MASS)=0.0;
758  }
759  }
760 
761  KRATOS_CATCH("")
762 }
763 
765 {
766  KRATOS_TRY
767 
768  ModelPart& r_model_part = BaseType::GetModelPart();
769  NodesArrayType& pNodes = r_model_part.Nodes();
770  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
771  #ifdef _OPENMP
772  int number_of_threads = omp_get_max_threads();
773  #else
774  int number_of_threads = 1;
775  #endif
776 
777  vector<unsigned int> node_partition;
778  CreatePartition(number_of_threads, pNodes.size(), node_partition);
779 
780  #pragma omp parallel for
781  for(int k=0; k<number_of_threads; k++)
782  {
783  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
784  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
785 
786  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
787  {
788  i->FastGetSolutionStepValue(PRESSUREAUX) = i->FastGetSolutionStepValue(NODAL_AREA)/i->FastGetSolutionStepValue(NODAL_MASS);
789  i->FastGetSolutionStepValue(EXTERNAL_PRESSURE) = i->FastGetSolutionStepValue(PRESSUREAUX) + i->FastGetSolutionStepValue(PRESSURE);
790  }
791  }
792 
793  KRATOS_CATCH("")
794 }
795 
796 
797 void SetToZeroMassAndArea(ProcessInfo& CurrentProcessInfo)
798 {
799  KRATOS_TRY
800 
801  ModelPart& r_model_part = BaseType::GetModelPart();
802  NodesArrayType& pNodes = r_model_part.Nodes();
803  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
804  #ifdef _OPENMP
805  int number_of_threads = omp_get_max_threads();
806  #else
807  int number_of_threads = 1;
808  #endif
809 
810  vector<unsigned int> node_partition;
811  CreatePartition(number_of_threads, pNodes.size(), node_partition);
812 
813  #pragma omp parallel for
814  for(int k=0; k<number_of_threads; k++)
815  {
816  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
817  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
818 
819  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
820  {
821  i->FastGetSolutionStepValue(NODAL_MASS) = 0.0;
822  i->FastGetSolutionStepValue(NODAL_AREA) = 0.0;
823  }
824  }
825 
826  KRATOS_CATCH("")
827 }
828 
829 void UpdateLoopForMassAndArea(ProcessInfo& CurrentProcessInfo)
830 {
831  KRATOS_TRY
832  KRATOS_CATCH("")
833 }
834 
835 
836 };
837 
838 
839  template<
840  class TSparseSpace,
841  class TDenseSpace,
842  class TLinearSolver>
843  class Fluid_Phase_PFEM2_Explicit_Strategy : public ImplicitSolvingStrategy<TSparseSpace,TDenseSpace,TLinearSolver>
844  {
845 
846  public:
847 
849 
851 
852  typedef typename BaseType::TDataType TDataType;
853 
854  typedef TSparseSpace SparseSpaceType;
855 
857 
859 
861 
863 
865 
867 
869 
871 
873 
875 
877 
879 
880  typedef ConditionsContainerType::iterator ConditionsContainerIterator;
881 
883 
885 
887 
888 
889  //typedef Element::Pointer ParticlePointer;
890  //typedef typename std::vector<ParticlePointer> ParticlePointerVector;
891  //typedef typename std::vector<ParticlePointer>::iterator ParticlePointerIterator;
892 
893  //typedef GlobalPointersVector<Element > ParticleWeakVector;
894  //typedef GlobalPointersVector<Element >::iterator ParticleWeakIterator;
895 
896 
897 
898 
901  const int dimension,
902  // const int damp_type,
903  //const double damping_ratio,
904  // const bool virtual_mass,
905  // const double contact_stiffness_ratio,
906  // const double max_delta_time,
907  //const bool CalculateReactions,
908  //const bool ComputeFemFemContact,
909  const bool MoveMeshFlag
910  //typename TLinearSolver::Pointer pNewLinearSolver,
911  //typename TSchemeType::Pointer pScheme,
912  //typename TBuilderAndSolverType::Pointer pNewBuilderAndSolver
913  )
914 
915  : ImplicitSolvingStrategy<TSparseSpace,TDenseSpace,TLinearSolver>(model_part, MoveMeshFlag)
916  {
917 
918  }
919 
921 
922 
923 
924  //********************************************
925  //********************************************
926  inline void CreatePartition(unsigned int number_of_threads, const int number_of_rows, vector<unsigned int>& partitions)
927  {
928  partitions.resize(number_of_threads+1);
929  int partition_size = number_of_rows / number_of_threads;
930  partitions[0] = 0;
931  partitions[number_of_threads] = number_of_rows;
932  for(unsigned int i = 1; i<number_of_threads; i++)
933  partitions[i] = partitions[i-1] + partition_size ;
934  }
935 
936 
938 {
939 
940  KRATOS_TRY
941  ModelPart& r_model_part = BaseType::GetModelPart();
942  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
943  ElementsArrayType& pElements = r_model_part.Elements();
944 
945 
946  typename ElementsArrayType::iterator it_begin = pElements.ptr_begin();
947  typename ElementsArrayType::iterator it_end = pElements.ptr_end();
948  for (ElementsArrayType::iterator it = it_begin; it != it_end; ++it)
949  {
950  it->AddExplicitContribution(CurrentProcessInfo);
951  }
952 
953 
954  KRATOS_CATCH("")
955 }
956 
957  //SPECIFIC FUNCTIONS FOR MY APPLICATION
958 void InitializeSolutionStep() override
959 {
960  KRATOS_TRY
961 
962  ModelPart& r_model_part = BaseType::GetModelPart();
963  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
964 
965  switch ( CurrentProcessInfo[FRACTIONAL_STEP] )
966  {
967  case 0:
968  {
970  break;
971  }
972  case 3:
973  {
975  break;
976  }
977 
978  default:
979  {
980  KRATOS_THROW_ERROR(std::logic_error,"Unexpected value for FRACTIONAL_STEP index: ", CurrentProcessInfo[FRACTIONAL_STEP]);
981  }
982  }
983  KRATOS_CATCH("")
984 }
985 
986 void FinalizeSolutionStep() override
987 {
988  KRATOS_TRY
989 
990  ModelPart& r_model_part = BaseType::GetModelPart();
991  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
992 
993  switch ( CurrentProcessInfo[FRACTIONAL_STEP] )
994  {
995  case 0:
996  {
998  break;
999  }
1000  case 3:
1001  {
1003  break;
1004  }
1005 
1006  default:
1007  {
1008  KRATOS_THROW_ERROR(std::logic_error,"Unexpected value for FRACTIONAL_STEP index: ", CurrentProcessInfo[FRACTIONAL_STEP]);
1009  }
1010  }
1011  KRATOS_CATCH("")
1012 }
1013 
1014 
1015 
1016 
1017 //VISCOUS ITERATIONS
1018 
1020 {
1021  KRATOS_TRY
1022 
1023  ModelPart& r_model_part = BaseType::GetModelPart();
1024  //ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
1025  NodesArrayType& pNodes = r_model_part.Nodes();
1026 
1027  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
1028  //const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
1029 
1030  #ifdef _OPENMP
1031  int number_of_threads = omp_get_max_threads();
1032  #else
1033  int number_of_threads = 1;
1034  #endif
1035 
1036  vector<unsigned int> node_partition;
1037  CreatePartition(number_of_threads, pNodes.size(), node_partition);
1038 
1039  #pragma omp parallel for
1040  for(int k=0; k<number_of_threads; k++)
1041  {
1042  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
1043  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
1044 
1045  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
1046  {
1047  noalias(i->FastGetSolutionStepValue(RHS)) = ZeroVector(3);
1048  i->FastGetSolutionStepValue(NODAL_MASS)=0.0;
1049  }
1050  }
1051 
1052  KRATOS_CATCH("")
1053 }
1054 
1056 {
1057  KRATOS_TRY
1058 
1059  ModelPart& r_model_part = BaseType::GetModelPart();
1060  //ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
1061  NodesArrayType& pNodes = r_model_part.Nodes();
1062 
1063  #ifdef _OPENMP
1064  int number_of_threads = omp_get_max_threads();
1065  #else
1066  int number_of_threads = 1;
1067  #endif
1068 
1069  vector<unsigned int> node_partition;
1070  CreatePartition(number_of_threads, pNodes.size(), node_partition);
1071 
1072  #pragma omp parallel for
1073  for(int k=0; k<number_of_threads; k++)
1074  {
1075  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
1076  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
1077 
1078  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
1079  {
1080  array_1d<double,3>& rhs = (i->FastGetSolutionStepValue(RHS));
1081  array_1d<double,3>& node_update_variable = (i->FastGetSolutionStepValue(WATER_VELOCITY));
1082  noalias(node_update_variable) = rhs/(i->FastGetSolutionStepValue(NODAL_MASS)) ;
1083  if (i->IsFixed(WATER_VELOCITY_X)==true)
1084  noalias(node_update_variable) = (i->FastGetSolutionStepValue(WATER_VELOCITY,1));
1085 
1086  }
1087  }
1088 
1089  KRATOS_CATCH("")
1090 }
1091 
1092 
1093 //PRESSURE ITERATIONS
1094 
1096 {
1097  KRATOS_TRY
1098 
1099  ModelPart& r_model_part = BaseType::GetModelPart();
1100  ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
1101  NodesArrayType& pNodes = r_model_part.Nodes();
1102 
1103  //const double delta_t = CurrentProcessInfo[DELTA_TIME];
1104  const int iteration_number = CurrentProcessInfo[NL_ITERATION_NUMBER];
1105 
1106  #ifdef _OPENMP
1107  int number_of_threads = omp_get_max_threads();
1108  #else
1109  int number_of_threads = 1;
1110  #endif
1111 
1112  vector<unsigned int> node_partition;
1113  CreatePartition(number_of_threads, pNodes.size(), node_partition);
1114 
1115  #pragma omp parallel for
1116  for(int k=0; k<number_of_threads; k++)
1117  {
1118  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
1119  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
1120 
1121  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
1122  {
1123 
1124  if (iteration_number == 1)
1125  noalias(i->FastGetSolutionStepValue(ACCELERATION)) = ZeroVector(3);
1126  else// if (iteration_number==2)
1127  noalias(i->FastGetSolutionStepValue(ACCELERATION)) = i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO); //second order cos we are in the second (or higher) iteration.
1128  //else if (iteration_number==3)
1129  // noalias(inode->FastGetSolutionStepValue(ACCELERATION)) = - inode->FastGetSolutionStepValue(ACCELERATION) + inode->FastGetSolutionStepValue(PRESS_PROJ_NO_RO);
1130 
1131  //noalias(i->FastGetSolutionStepValue(PRESS_PROJ)) = ZeroVector(3);
1132  noalias(i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO)) = ZeroVector(3);
1133 
1134 
1135  //noalias(in->GetSolutionStepValue(PRESSURE,1))=in->FastGetSolutionStepValue(PRESSURE);
1136  //noalias(in->FastGetSolutionStepValue(PRESSURE)) = 0.0;
1137  //i->FastGetSolutionStepValue(NODAL_AREA)=0.0;
1138  i->FastGetSolutionStepValue(NODAL_MASS)=0.0;
1139  //i->FastGetSolutionStepValue(MASS)=ZeroVector(3);
1140  }
1141  }
1142 
1143  KRATOS_CATCH("")
1144 }
1145 
1147 {
1148  KRATOS_TRY
1149 
1150  ModelPart& r_model_part = BaseType::GetModelPart();
1151  //ProcessInfo& CurrentProcessInfo = r_model_part.GetProcessInfo();
1152  NodesArrayType& pNodes = r_model_part.Nodes();
1153  //const double factor = CurrentProcessInfo.GetValue(DELTA_TIME); //included in factor
1154  #ifdef _OPENMP
1155  int number_of_threads = omp_get_max_threads();
1156  #else
1157  int number_of_threads = 1;
1158  #endif
1159 
1160  vector<unsigned int> node_partition;
1161  CreatePartition(number_of_threads, pNodes.size(), node_partition);
1162 
1163  #pragma omp parallel for
1164  for(int k=0; k<number_of_threads; k++)
1165  {
1166  typename NodesArrayType::iterator i_begin=pNodes.ptr_begin()+node_partition[k];
1167  typename NodesArrayType::iterator i_end=pNodes.ptr_begin()+node_partition[k+1];
1168 
1169  for(ModelPart::NodeIterator i=i_begin; i!= i_end; ++i)
1170  {
1171  //normalizing variables:
1172  array_1d<double,3>& press_proj_no_ro = (i->FastGetSolutionStepValue(PRESS_PROJ_NO_RO));
1173  //array_1d<double,3>& press_proj_stabilization = (i->FastGetSolutionStepValue(PRESS_PROJ));
1174  double& mass = (i->FastGetSolutionStepValue(NODAL_MASS));
1175  //array_1d<double,3>& vectorial_mass=(i->FastGetSolutionStepValue(MASS));
1176  //double& area = (i->FastGetSolutionStepValue(NODAL_AREA));
1177  press_proj_no_ro /= mass;
1178  //press_proj_no_ro(0) /= vectorial_mass(0);
1179  //press_proj_no_ro(1) /= vectorial_mass(1);
1180  //press_proj_no_ro(2) /= vectorial_mass(2)+1e-20;
1181  //press_proj_stabilization /= area;
1182 
1183  //updating acceleration
1184  array_1d<double,3>& acceleration = (i->FastGetSolutionStepValue(ACCELERATION));
1185  noalias(acceleration) -= (press_proj_no_ro);
1186  //updating variable
1187  array_1d<double,3>& velocity = (i->FastGetSolutionStepValue(WATER_VELOCITY));
1188  noalias(velocity) += (acceleration) ;
1189 
1190  i->FastGetSolutionStepValue(PREVIOUS_ITERATION_PRESSURE)=i->FastGetSolutionStepValue(WATER_PRESSURE);
1191  //i->FastGetSolutionStepValue(FRACT_VEL)=i->FastGetSolutionStepValue(VELOCITY);
1192  if (i->IsFixed(WATER_VELOCITY_X)==true)
1193  noalias(velocity) = (i->FastGetSolutionStepValue(WATER_VELOCITY,1));
1194  }
1195  }
1196 
1197  KRATOS_CATCH("")
1198 }
1199 
1200 
1201 };
1202 
1203 } /* namespace Kratos.*/
1204 #endif /* KRATOS_RESIDUALBASED_CENTRAL_DIFERENCES_STRATEGY */
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Definition: explicit_strategy.h:844
BaseType::TSystemMatrixType TSystemMatrixType
Definition: explicit_strategy.h:864
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: explicit_strategy.h:876
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: explicit_strategy.h:870
BaseType::TDataType TDataType
Definition: explicit_strategy.h:852
BaseType::TBuilderAndSolverType TBuilderAndSolverType
Definition: explicit_strategy.h:856
BaseType::DofsArrayType DofsArrayType
Definition: explicit_strategy.h:860
BaseType::TSystemVectorType TSystemVectorType
Definition: explicit_strategy.h:866
Fluid_Phase_PFEM2_Explicit_Strategy(ModelPart &model_part, const int dimension, const bool MoveMeshFlag)
Definition: explicit_strategy.h:899
ModelPart::NodesContainerType NodesArrayType
Definition: explicit_strategy.h:872
void AssembleLoop()
Definition: explicit_strategy.h:937
void CreatePartition(unsigned int number_of_threads, const int number_of_rows, vector< unsigned int > &partitions)
Definition: explicit_strategy.h:926
BaseType::TSystemMatrixPointerType TSystemMatrixPointerType
Definition: explicit_strategy.h:882
ModelPart::ElementsContainerType ElementsArrayType
Definition: explicit_strategy.h:874
BaseType::TSchemeType TSchemeType
Definition: explicit_strategy.h:858
TSparseSpace SparseSpaceType
Definition: explicit_strategy.h:854
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: explicit_strategy.h:850
BaseType::TSystemVectorPointerType TSystemVectorPointerType
Definition: explicit_strategy.h:884
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: explicit_strategy.h:868
virtual ~Fluid_Phase_PFEM2_Explicit_Strategy()
Definition: explicit_strategy.h:920
void UpdateLoopForViscousIterationsWithNormalization()
Definition: explicit_strategy.h:1055
void UpdateLoopForPressureIterationsWithNormalization()
Definition: explicit_strategy.h:1146
ModelPart::PropertiesType PropertiesType
Definition: explicit_strategy.h:886
KRATOS_CLASS_POINTER_DEFINITION(Fluid_Phase_PFEM2_Explicit_Strategy)
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: explicit_strategy.h:986
Element::DofsVectorType DofsVectorType
Definition: explicit_strategy.h:862
void SetToZeroVariablesInViscousIterations()
Definition: explicit_strategy.h:1019
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: explicit_strategy.h:958
ModelPart::ConditionsContainerType::ContainerType ConditionsContainerType
Definition: explicit_strategy.h:878
ConditionsContainerType::iterator ConditionsContainerIterator
Definition: explicit_strategy.h:880
void SetToZeroVariablesInPresureIterations()
Definition: explicit_strategy.h:1095
Implicit solving strategy base class This is the base class from which we will derive all the implici...
Definition: implicit_solving_strategy.h:61
Scheme< TSparseSpace, TDenseSpace > TSchemeType
Definition: implicit_solving_strategy.h:82
BaseType::NodesArrayType NodesArrayType
Definition: implicit_solving_strategy.h:92
BuilderAndSolver< TSparseSpace, TDenseSpace, TLinearSolver > TBuilderAndSolverType
Definition: implicit_solving_strategy.h:84
BaseType::ElementsArrayType ElementsArrayType
Definition: implicit_solving_strategy.h:94
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ConditionsContainerType ConditionsContainerType
Condintions container. A vector set of Conditions with their Id's as key.
Definition: model_part.h:183
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
Definition: explicit_strategy.h:96
ModelPart::ConditionsContainerType::ContainerType ConditionsContainerType
Definition: explicit_strategy.h:130
void NormalizePressureProjection(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:614
BaseType::DofsArrayType DofsArrayType
Definition: explicit_strategy.h:112
ConditionsContainerType::iterator ConditionsContainerIterator
Definition: explicit_strategy.h:132
ImplicitSolvingStrategy< TSparseSpace, TDenseSpace, TLinearSolver > BaseType
Definition: explicit_strategy.h:102
BaseType::TSystemMatrixType TSystemMatrixType
Definition: explicit_strategy.h:116
ModelPart::ConditionsContainerType ConditionsArrayType
Definition: explicit_strategy.h:128
void SetToZeroVariablesForPressure(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:489
PFEM2_Explicit_Strategy(ModelPart &model_part, const int dimension, const bool MoveMeshFlag)
Definition: explicit_strategy.h:151
BaseType::TSystemVectorPointerType TSystemVectorPointerType
Definition: explicit_strategy.h:136
void UpdateLoopForViscousIterationsWithNormalization(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:358
void InitializeSolutionStep() override
Performs all the required operations that should be done (for each step) before solving the solution ...
Definition: explicit_strategy.h:210
void SetToZeroVariablesForVolumetricStrain(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:454
ModelPart::ElementsContainerType ElementsArrayType
Definition: explicit_strategy.h:126
ModelPart::NodesContainerType NodesArrayType
Definition: explicit_strategy.h:124
void AssembleLoop()
Definition: explicit_strategy.h:189
virtual ~PFEM2_Explicit_Strategy()
Definition: explicit_strategy.h:172
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: explicit_strategy.h:120
BaseType::TDataType TDataType
Definition: explicit_strategy.h:104
BaseType::TSystemVectorType TSystemVectorType
Definition: explicit_strategy.h:118
ModelPart::PropertiesType PropertiesType
Definition: explicit_strategy.h:138
void UpdateLoopForVolumetricStrain(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:660
void FinalizeSolutionStep() override
Performs all the required operations that should be done (for each step) after solving the solution s...
Definition: explicit_strategy.h:265
void UpdateLoopForPressure(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:692
BaseType::TSystemMatrixPointerType TSystemMatrixPointerType
Definition: explicit_strategy.h:134
Element::DofsVectorType DofsVectorType
Definition: explicit_strategy.h:114
void UpdateLoopForPressureViscousCorrection(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:764
void SetToZeroVariablesInPresureViscousCorrection(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:729
void CreatePartition(unsigned int number_of_threads, const int number_of_rows, vector< unsigned int > &partitions)
Definition: explicit_strategy.h:178
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: explicit_strategy.h:122
void UpdateLoopForMassAndArea(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:829
KRATOS_CLASS_POINTER_DEFINITION(PFEM2_Explicit_Strategy)
BaseType::TSchemeType TSchemeType
Definition: explicit_strategy.h:110
void SetToZeroVariablesInPresureProjection(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:575
TSparseSpace SparseSpaceType
Definition: explicit_strategy.h:106
void SetToZeroMassAndArea(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:797
void SetToZeroVariablesInPresureIterations(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:398
BaseType::TBuilderAndSolverType TBuilderAndSolverType
Definition: explicit_strategy.h:108
void SetToZeroVariablesInViscousIterations(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:323
void UpdateLoopForPressureIterationsWithNormalization(ProcessInfo &CurrentProcessInfo)
Definition: explicit_strategy.h:524
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
TSparseSpace::VectorPointerType TSystemVectorPointerType
Definition: solving_strategy.h:77
TDenseSpace::MatrixType LocalSystemMatrixType
Definition: solving_strategy.h:79
TSparseSpace::DataType TDataType
Definition: solving_strategy.h:69
ModelPart & GetModelPart()
Operations to get the pointer to the model.
Definition: solving_strategy.h:350
TSparseSpace::MatrixType TSystemMatrixType
Definition: solving_strategy.h:71
TSparseSpace::MatrixPointerType TSystemMatrixPointerType
Definition: solving_strategy.h:75
bool MoveMeshFlag()
This function returns the flag that says if the mesh is moved.
Definition: solving_strategy.h:290
TSparseSpace::VectorType TSystemVectorType
Definition: solving_strategy.h:73
TDenseSpace::VectorType LocalSystemVectorType
Definition: solving_strategy.h:81
#define KRATOS_THROW_ERROR(ExceptionType, ErrorMessage, MoreInfo)
Definition: define.h:77
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
float velocity
Definition: PecletTest.py:54
model_part
Definition: face_heat.py:14
tuple acceleration
Definition: generate_droplet_dynamics.py:117
rhs
Definition: generate_frictional_mortar_condition.py:297
int dimension
Definition: isotropic_damage_automatic_differentiation.py:123
int k
Definition: quadrature.py:595
integer i
Definition: TensorModule.f:17
Configure::ContainerType ContainerType
Definition: transfer_utility.h:247