595 assert(mProblemDimension==DIM+1);
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;
609 iter != mrQuadMesh.GetElementIteratorEnd();
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] );
623 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
624 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
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] );
631 left_val = mCurrentSolution[mProblemDimension*vertex_0_global_index + DIM];
632 right_val = mCurrentSolution[mProblemDimension*vertex_1_global_index + DIM];
636 mCurrentSolution[mProblemDimension*global_index + DIM] = 0.5 * (left_val + right_val);
674 std::vector<unsigned> rows;
675 std::vector<double> values;
678 bool applySymmetrically = (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY) && symmetricProblem;
680 if (applySymmetrically)
682 if (mDirichletBoundaryConditionsVector ==
nullptr)
684 VecDuplicate(mResidualVector, &mDirichletBoundaryConditionsVector);
695 for (
unsigned i=0; i<mrProblemDefinition.rGetDirichletNodes().size(); i++)
697 unsigned node_index = mrProblemDefinition.rGetDirichletNodes()[i];
699 for (
unsigned j=0; j<DIM; j++)
701 double dirichlet_val = mrProblemDefinition.rGetDirichletNodeValues()[i](j);
706 unsigned dof_index = mProblemDimension*node_index+j;
708 if (type == LINEAR_PROBLEM)
719 val = mCurrentSolution[dof_index] - dirichlet_val;
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];
759 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
764 if (applySymmetrically)
781 if (type!=LINEAR_PROBLEM)
783 for (
unsigned i=0; i<rows.size(); i++)
789 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
791 for (
unsigned i=0; i<rows.size(); i++)
801 assert(mCompressibilityType==INCOMPRESSIBLE);
804 VecGetOwnershipRange(mResidualVector, &lo, &hi);
806 for (
unsigned i=0; i<mrQuadMesh.GetNumNodes(); i++)
808 if (mrQuadMesh.GetNode(i)->IsInternal())
810 unsigned row = (DIM+1)*i + DIM;
811 if (lo <= (
int)row && (int)row < hi)
813 if (type!=LINEAR_PROBLEM)
817 if (type!=NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY)
819 double rhs_vector_val = type==LINEAR_PROBLEM ? 0.0 : mCurrentSolution[row]-0.0;
833 Vec template_vec = mrQuadMesh.GetDistributedVectorFactory()->CreateVec(mProblemDimension);
838 VecDuplicate(template_vec, &mResidualVector);
839 VecDuplicate(mResidualVector, &mLinearSystemRhsVector);
842 mDirichletBoundaryConditionsVector =
nullptr;
850 VecGetOwnershipRange(mResidualVector, &lo, &hi);
857 unsigned num_non_zeros = std::min(75u, mNumDofs);
859 PetscTools::SetupMat(mSystemLhsMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
860 PetscTools::SetupMat(mPreconditionMatrix, mNumDofs, mNumDofs, num_non_zeros, local_size, local_size);
869 int* num_non_zeros_each_row =
new int[mNumDofs];
870 for (
unsigned i=0; i<mNumDofs; i++)
872 num_non_zeros_each_row[i] = 0;
876 iter != mrQuadMesh.GetNodeIteratorEnd();
882 unsigned num_non_zeros_upper_bound = 4 + 30*iter->GetNumContainingElements();
884 num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, mNumDofs);
886 unsigned i = iter->GetIndex();
888 num_non_zeros_each_row[mProblemDimension*i + 0] = num_non_zeros_upper_bound;
889 num_non_zeros_each_row[mProblemDimension*i + 1] = num_non_zeros_upper_bound;
890 num_non_zeros_each_row[mProblemDimension*i + 2] = num_non_zeros_upper_bound;
892 if (mCompressibilityType==INCOMPRESSIBLE)
894 if (!iter->IsInternal())
896 num_non_zeros_each_row[mProblemDimension*i + 3] = num_non_zeros_upper_bound;
900 num_non_zeros_each_row[mProblemDimension*i + 3] = 1;
919#if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2)
920 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mSystemLhsMatrix);
921 MatCreate(PETSC_COMM_WORLD,local_size,local_size,mNumDofs,mNumDofs,&mPreconditionMatrix);
923 MatCreate(PETSC_COMM_WORLD,&mSystemLhsMatrix);
924 MatCreate(PETSC_COMM_WORLD,&mPreconditionMatrix);
925 MatSetSizes(mSystemLhsMatrix,local_size,local_size,mNumDofs,mNumDofs);
926 MatSetSizes(mPreconditionMatrix,local_size,local_size,mNumDofs,mNumDofs);
931 MatSetType(mSystemLhsMatrix, MATSEQAIJ);
932 MatSetType(mPreconditionMatrix, MATSEQAIJ);
933 MatSeqAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
934 MatSeqAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row);
938 int* num_non_zeros_each_row_in_diag =
new int[local_size];
939 int* num_non_zeros_each_row_off_diag =
new int[local_size];
940 for (
unsigned i=0; i<
unsigned(local_size); i++)
942 num_non_zeros_each_row_in_diag[i] = num_non_zeros_each_row[lo+i];
943 num_non_zeros_each_row_off_diag[i] = num_non_zeros_each_row[lo+i];
945 if (num_non_zeros_each_row_in_diag[i] > local_size)
947 num_non_zeros_each_row_in_diag[i] = local_size;
951 MatSetType(mSystemLhsMatrix, MATMPIAIJ);
952 MatSetType(mPreconditionMatrix, MATMPIAIJ);
953 MatMPIAIJSetPreallocation(mSystemLhsMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
954 MatMPIAIJSetPreallocation(mPreconditionMatrix, PETSC_DEFAULT, num_non_zeros_each_row_in_diag, PETSC_DEFAULT, num_non_zeros_each_row_off_diag);
957 MatSetFromOptions(mSystemLhsMatrix);
958 MatSetFromOptions(mPreconditionMatrix);
959#if (PETSC_VERSION_MAJOR == 3)
960 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
961 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
963 MatSetOption(mSystemLhsMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
964 MatSetOption(mPreconditionMatrix, MAT_IGNORE_OFF_PROC_ENTRIES);
974 delete [] num_non_zeros_each_row;