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");
377 template<
unsigned DIM>
380 std::string outputDirectory,
381 CompressibilityType compressibilityType)
382 : mrQuadMesh(rQuadMesh),
383 mrProblemDefinition(rProblemDefinition),
384 mOutputDirectory(outputDirectory),
385 mpOutputFileHandler(NULL),
386 mpQuadratureRule(NULL),
387 mpBoundaryQuadratureRule(NULL),
388 mCompressibilityType(compressibilityType),
389 mResidualVector(NULL),
390 mSystemLhsMatrix(NULL),
391 mPreconditionMatrix(NULL)
393 assert(DIM==2 || DIM==3);
399 if(p_quad_mesh == NULL && p_distributed_quad_mesh == NULL)
401 EXCEPTION(
"Continuum mechanics solvers require a quadratic mesh");
432 template<
unsigned DIM>
435 if (mpOutputFileHandler)
437 delete mpOutputFileHandler;
440 if (mpQuadratureRule)
442 delete mpQuadratureRule;
443 delete mpBoundaryQuadratureRule;
454 if (mDirichletBoundaryConditionsVector)
460 template<
unsigned DIM>
462 std::string fileExtension,
473 std::stringstream file_name;
475 file_name << fileName;
476 if (counterToAppend >= 0)
478 file_name <<
"_" << counterToAppend;
480 file_name <<
"." << fileExtension;
482 out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
484 std::vector<c_vector<double,DIM> >& r_spatial_solution = rGetSpatialSolution();
485 for (
unsigned i=0; i<r_spatial_solution.size(); i++)
491 for (
unsigned j=0; j<DIM; j++)
493 *p_file << r_spatial_solution[i](j) <<
" ";
502 template<
unsigned DIM>
513 std::stringstream file_name;
515 file_name <<
"pressure";
516 if (counterToAppend >= 0)
518 file_name <<
"_" << counterToAppend;
522 out_stream p_file = mpOutputFileHandler->OpenOutputFile(file_name.str());
524 std::vector<double>& r_pressure = rGetPressures();
525 for (
unsigned i=0; i<r_pressure.size(); i++)
527 for (
unsigned j=0; j<DIM; j++)
529 *p_file << mrQuadMesh.GetNode(i)->rGetLocation()[j] <<
" ";
532 *p_file << r_pressure[i] <<
"\n";
540 template<
unsigned DIM>
543 if (writeOutput && (mOutputDirectory==
""))
545 EXCEPTION(
"Can't write output if no output directory was given in constructor");
547 mWriteOutput = writeOutput;
551 template<
unsigned DIM>
554 if (this->mOutputDirectory==
"")
556 EXCEPTION(
"No output directory was given so no output was written, cannot convert to VTK format");
561 mesh_writer.
AddPointData(spatialSolutionName, this->rGetSpatialSolution());
563 if (mCompressibilityType==INCOMPRESSIBLE)
569 std::vector<double> element_attribute;
571 iter != this->mrQuadMesh.GetElementIteratorEnd();
574 element_attribute.push_back(iter->GetAttribute());
576 mesh_writer.
AddCellData(
"Attribute", element_attribute);
583 template<
unsigned DIM>
586 assert(mProblemDimension==DIM+1);
588 mPressureSolution.clear();
589 mPressureSolution.resize(mrQuadMesh.GetNumNodes());
591 for (
unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
593 mPressureSolution[i] = mCurrentSolution[mProblemDimension*i + DIM];
595 return mPressureSolution;
599 template<
unsigned DIM>
602 assert(mProblemDimension==DIM+1);
605 unsigned internal_nodes_2d[3] = {3,4,5};
606 unsigned neighbouring_vertices_2d[3][2] = { {1,2}, {2,0}, {0,1} };
609 unsigned internal_nodes_3d[6] = {4,5,6,7,8,9};
610 unsigned neighbouring_vertices_3d[6][2] = { {0,1}, {1,2}, {0,2}, {0,3}, {1,3}, {2,3} };
612 unsigned num_internal_nodes_per_element = DIM==2 ? 3 : 6;
616 iter != mrQuadMesh.GetElementIteratorEnd();
619 for(
unsigned i=0; i<num_internal_nodes_per_element; i++)
621 unsigned global_index;
627 global_index = iter->GetNodeGlobalIndex( internal_nodes_2d[i] );
628 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][0] );
629 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_2d[i][1] );
630 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
631 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
635 global_index = iter->GetNodeGlobalIndex( internal_nodes_3d[i] );
636 unsigned vertex_0_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][0] );
637 unsigned vertex_1_global_index =iter->GetNodeGlobalIndex( neighbouring_vertices_3d[i][1] );
638 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
639 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
643 mCurrentSolution[mProblemDimension*global_index + DIM] = 0.5 * (left_val + right_val);
678 template<
unsigned DIM>
681 std::vector<unsigned> rows;
682 std::vector<double> values;
685 bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
687 if (applySymmetrically)
689 if(mDirichletBoundaryConditionsVector==NULL)
691 VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
702 for (
unsigned i=0; i<mrProblemDefinition.rGetDirichletNodes().size(); i++)
704 unsigned node_index = mrProblemDefinition.rGetDirichletNodes()[i];
706 for (
unsigned j=0; j<DIM; j++)
708 double dirichlet_val = mrProblemDefinition.rGetDirichletNodeValues()[i](j);
713 unsigned dof_index = mProblemDimension*node_index+j;
715 if(type == LINEAR_PROBLEM)
726 val = mCurrentSolution[dof_index] - dirichlet_val;
728 rows.push_back(dof_index);
729 values.push_back(val);
738 if (applySymmetrically)
741 for (
unsigned i=0; i<rows.size(); i++)
743 unsigned col = rows[i];
744 double minus_value = -values[i];
765 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
770 if (applySymmetrically)
785 if (type!=LINEAR_PROBLEM)
787 for (
unsigned i=0; i<rows.size(); i++)
793 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
795 for (
unsigned i=0; i<rows.size(); i++)
803 template<
unsigned DIM>
806 assert(mCompressibilityType==INCOMPRESSIBLE);
809 VecGetOwnershipRange(mResidualVector, &lo, &hi);
811 for(
unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
813 if(mrQuadMesh.GetNode(i)->IsInternal())
815 unsigned row = (DIM+1)*i + DIM;
816 if(lo <= (
int)row && (int)row < hi)
818 if(type!=LINEAR_PROBLEM)
822 if(type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
824 double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 : mCurrentSolution[row]-0.0;
837 template<
unsigned DIM>
840 Vec template_vec = mrQuadMesh.GetDistributedVectorFactory()->CreateVec(mProblemDimension);
845 VecDuplicate(template_vec, &mResidualVector);
846 VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
849 mDirichletBoundaryConditionsVector = NULL;
857 VecGetOwnershipRange(mResidualVector, &lo, &hi);
864 unsigned num_non_zeros = std::min(75u, mNumDofs);
866 PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
867 PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
876 int* num_non_zeros_each_row =
new int[mNumDofs];
877 for (
unsigned i=0; i<mNumDofs; i++)
879 num_non_zeros_each_row[i] = 0;
883 iter != mrQuadMesh.GetNodeIteratorEnd();
889 unsigned num_non_zeros_upper_bound = 4 + 30*iter->GetNumContainingElements();
891 num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
893 unsigned i = iter->GetIndex();
895 num_non_zeros_each_row[mProblemDimension*i + 0] = num_non_zeros_upper_bound;
896 num_non_zeros_each_row[mProblemDimension*i + 1] = num_non_zeros_upper_bound;
897 num_non_zeros_each_row[mProblemDimension*i + 2] = num_non_zeros_upper_bound;
899 if (mCompressibilityType==INCOMPRESSIBLE)
901 if(!iter->IsInternal())
903 num_non_zeros_each_row[mProblemDimension*i + 3] = num_non_zeros_upper_bound;
907 num_non_zeros_each_row[mProblemDimension*i + 3] = 1;
926 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
927 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mSystemLhsMatrix);
928 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mPreconditionMatrix);
930 MatCreate(PETSC_COMM_WORLD,&mSystemLhsMatrix);
931 MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
932 MatSetSizes(mSystemLhsMatrix,local_size,local_size,mNumDofs,mNumDofs);
933 MatSetSizes(mPreconditionMatrix,local_size,local_size,mNumDofs,mNumDofs);
938 MatSetType(mSystemLhsMatrix, MATSEQAIJ);
939 MatSetType(mPreconditionMatrix, MATSEQAIJ);
940 MatSeqAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
941 MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
945 int* num_non_zeros_each_row_in_diag =
new int[local_size];
946 int* num_non_zeros_each_row_off_diag =
new int[local_size];
947 for (
unsigned i=0; i<
unsigned(local_size); i++)
949 num_non_zeros_each_row_in_diag[i] = num_non_zeros_each_row[lo+i];
950 num_non_zeros_each_row_off_diag[i] = num_non_zeros_each_row[lo+i];
952 if(num_non_zeros_each_row_in_diag[i] > local_size)
954 num_non_zeros_each_row_in_diag[i] = local_size;
958 MatSetType(mSystemLhsMatrix, MATMPIAIJ);
959 MatSetType(mPreconditionMatrix, MATMPIAIJ);
960 MatMPIAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
961 MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
964 MatSetFromOptions(mSystemLhsMatrix);
965 MatSetFromOptions(mPreconditionMatrix);
966 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
967 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
968 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
970 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
971 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
981 delete [] num_non_zeros_each_row;
984 #endif // ABSTRACTCONTINUUMMECHANICSSOLVER_HPP_
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
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
void AllocateMatrixMemory()
void WriteCurrentSpatialSolution(std::string fileName, std::string fileExtension, int counterToAppend=-1)
#define EXCEPTION(message)
virtual ~AbstractContinuumMechanicsSolver()
std::vector< double > & rGetPressures()
Vec mDirichletBoundaryConditionsVector
void AddCellData(std::string name, std::vector< double > data)
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)
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)