36 #ifndef ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_ 37 #define ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_ 39 #include "ContinuumMechanicsProblemDefinition.hpp" 40 #include "CompressibilityType.hpp" 41 #include "LinearSystem.hpp" 42 #include "OutputFileHandler.hpp" 43 #include "ReplicatableVector.hpp" 44 #include "AbstractTetrahedralMesh.hpp" 45 #include "QuadraticMesh.hpp" 46 #include "DistributedQuadraticMesh.hpp" 47 #include "Warnings.hpp" 48 #include "PetscException.hpp" 49 #include "GaussianQuadratureRule.hpp" 51 #include "MechanicsEventHandler.hpp" 52 #include "CommandLineArguments.hpp" 53 #include "VtkMeshWriter.hpp" 61 typedef enum _ApplyDirichletBcsType
64 NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY,
65 NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING
66 } ApplyDirichletBcsType;
86 template<
unsigned DIM>
300 std::string outputDirectory,
301 CompressibilityType compressibilityType);
347 void CreateVtkOutput(std::string spatialSolutionName=
"Spatial solution");
375 template<
unsigned DIM>
378 std::string outputDirectory,
379 CompressibilityType compressibilityType)
391 assert(DIM==2 || DIM==3);
397 if ((p_quad_mesh ==
nullptr) && (p_distributed_quad_mesh ==
nullptr))
399 EXCEPTION(
"Continuum mechanics solvers require a quadratic mesh");
429 template<
unsigned DIM>
457 template<
unsigned DIM>
459 std::string fileExtension,
470 std::stringstream file_name;
472 file_name << fileName;
473 if (counterToAppend >= 0)
475 file_name <<
"_" << counterToAppend;
477 file_name <<
"." << fileExtension;
482 for (
unsigned i=0; i<r_spatial_solution.size(); i++)
488 for (
unsigned j=0; j<DIM; j++)
490 *p_file << r_spatial_solution[i](j) <<
" ";
499 template<
unsigned DIM>
508 std::stringstream file_name;
510 file_name <<
"pressure";
511 if (counterToAppend >= 0)
513 file_name <<
"_" << counterToAppend;
520 for (
unsigned i = 0; i < r_pressure.size(); i++)
522 for (
unsigned j = 0; j < DIM; j++)
527 *p_file << r_pressure[i] <<
"\n";
535 template<
unsigned DIM>
540 EXCEPTION(
"Can't write output if no output directory was given in constructor");
546 template<
unsigned DIM>
551 EXCEPTION(
"No output directory was given so no output was written, cannot convert to VTK format");
564 std::vector<double> element_attribute;
569 element_attribute.push_back(iter->GetAttribute());
571 mesh_writer.
AddCellData(
"Attribute", element_attribute);
577 template<
unsigned DIM>
592 template<
unsigned DIM>
598 unsigned internal_nodes_2d[3] = {3,4,5};
599 unsigned neighbouring_vertices_2d[3][2] = { {1,2}, {2,0}, {0,1} };
602 unsigned internal_nodes_3d[6] = {4,5,6,7,8,9};
603 unsigned neighbouring_vertices_3d[6][2] = { {0,1}, {1,2}, {0,2}, {0,3}, {1,3}, {2,3} };
605 unsigned num_internal_nodes_per_element = DIM==2 ? 3 : 6;
612 for (
unsigned i=0; i<num_internal_nodes_per_element; i++)
614 unsigned global_index;
620 global_index = iter->GetNodeGlobalIndex( internal_nodes_2d[i] );
621 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][0] );
622 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][1] );
628 global_index = iter->GetNodeGlobalIndex( internal_nodes_3d[i] );
629 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][0] );
630 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][1] );
671 template<
unsigned DIM>
674 std::vector<unsigned> rows;
675 std::vector<double> values;
678 bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
680 if (applySymmetrically)
699 for (
unsigned j=0; j<DIM; j++)
708 if (type == LINEAR_PROBLEM)
721 rows.push_back(dof_index);
722 values.push_back(val);
731 if (applySymmetrically)
734 for (
unsigned i=0; i<rows.size(); i++)
736 unsigned col = rows[i];
737 double minus_value = -values[i];
758 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
763 if (applySymmetrically)
778 if (type!=LINEAR_PROBLEM)
780 for (
unsigned i=0; i<rows.size(); i++)
786 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
788 for (
unsigned i=0; i<rows.size(); i++)
795 template<
unsigned DIM>
807 unsigned row = (DIM+1)*i + DIM;
808 if (lo <= (
int)row && (int)row < hi)
810 if (type!=LINEAR_PROBLEM)
814 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
816 double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 :
mCurrentSolution[row]-0.0;
827 template<
unsigned DIM>
854 unsigned num_non_zeros = std::min(75u,
mNumDofs);
866 int* num_non_zeros_each_row =
new int[
mNumDofs];
869 num_non_zeros_each_row[i] = 0;
879 unsigned num_non_zeros_upper_bound = 4 + 30*iter->GetNumContainingElements();
881 num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
883 unsigned i = iter->GetIndex();
891 if (!iter->IsInternal())
916 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 917 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&
mSystemLhsMatrix);
930 MatSeqAIJSetPreallocation(
mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
935 int* num_non_zeros_each_row_in_diag =
new int[local_size];
936 int* num_non_zeros_each_row_off_diag =
new int[local_size];
937 for (
unsigned i=0; i<
unsigned(local_size); i++)
939 num_non_zeros_each_row_in_diag[i] = num_non_zeros_each_row[lo+i];
940 num_non_zeros_each_row_off_diag[i] = num_non_zeros_each_row[lo+i];
942 if (num_non_zeros_each_row_in_diag[i] > local_size)
944 num_non_zeros_each_row_in_diag[i] = local_size;
950 MatMPIAIJSetPreallocation(
mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
951 MatMPIAIJSetPreallocation(
mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
956 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x 971 delete [] num_non_zeros_each_row;
974 #endif // ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_ ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
virtual DistributedVectorFactory * GetDistributedVectorFactory()
void ApplyDirichletBoundaryConditions(ApplyDirichletBcsType type, bool symmetricProblem)
GaussianQuadratureRule< DIM > * mpQuadratureRule
void WriteCurrentPressureSolution(int counterToAppend=-1)
Vec mLinearSystemRhsVector
AbstractTetrahedralMesh< DIM, DIM > & mrQuadMesh
virtual std::vector< c_vector< double, DIM > > & rGetSpatialSolution()=0
ElementIterator GetElementIteratorEnd()
void AllocateMatrixMemory()
void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1)
Node< SPACE_DIM > * GetNode(unsigned index) const
#define EXCEPTION(message)
virtual ~AbstractContinuumMechanicsSolver()
std::vector< double > & rGetPressures()
Vec mDirichletBoundaryConditionsVector
void AddCellData(std::string name, std::vector< double > data)
NodeIterator GetNodeIteratorEnd()
std::vector< c_vector< double, DIM > > mSpatialSolution
virtual unsigned GetNumNodes() const
AbstractContinuumMechanicsSolver(AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, ContinuumMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory, CompressibilityType compressibilityType)
void AddIdentityBlockForDummyPressureVariables(ApplyDirichletBcsType type)
bool OptionExists(const std::string &rOption)
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
out_stream OpenOutputFile(const std::string &rFileName, std::ios_base::openmode mode=std::ios::out|std::ios::trunc) const
void RemovePressureDummyValuesThroughLinearInterpolation()
std::vector< double > & rGetCurrentSolution()
ContinuumMechanicsProblemDefinition< DIM > & mrProblemDefinition
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
OutputFileHandler * mpOutputFileHandler
void CreateVtkOutput(std::string spatialSolutionName="Spatial solution")
CompressibilityType mCompressibilityType
std::vector< double > mCurrentSolution
std::vector< double > mPressureSolution
unsigned mProblemDimension
static CommandLineArguments * Instance()
GaussianQuadratureRule< DIM-1 > * mpBoundaryQuadratureRule
std::string mOutputDirectory
void AddPointData(std::string name, std::vector< double > data)
void SetWriteOutput(bool writeOutput=true)