void PetscVecTools::SetupInterleavedVectorScatterGather(Vec interleavedVec, VecScatter& rFirstVariableScatterContext, VecScatter& rSecondVariableScatterContext) { PetscInt num_rows, num_local_rows; VecGetSize(interleavedVec, &num_rows); VecGetLocalSize(interleavedVec, &num_local_rows); IS A11_rows, A22_rows; ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 2, &A11_rows); ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 1, 2, &A22_rows); IS all_vector; ISCreateStride(PETSC_COMM_WORLD, num_rows/2, 0, 1, &all_vector); unsigned subvector_num_rows = num_rows/2; unsigned subvector_local_rows = num_local_rows/2; Vec x1_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows); Vec x2_subvector = PetscTools::CreateVec(subvector_num_rows, subvector_local_rows); VecScatterCreate(interleavedVec, A11_rows, x1_subvector, all_vector, &rFirstVariableScatterContext); VecScatterCreate(interleavedVec, A22_rows, x2_subvector, all_vector, &rSecondVariableScatterContext); PetscTools::Destroy(x1_subvector); PetscTools::Destroy(x2_subvector); ISDestroy(PETSC_DESTROY_PARAM(A11_rows)); ISDestroy(PETSC_DESTROY_PARAM(A22_rows)); ISDestroy(PETSC_DESTROY_PARAM(all_vector)); }
void TestInterleavedVecScatter() { // Vectors will be twice PROBLEM_SIZE, since this is to be used in bidomain code. const unsigned PROBLEM_SIZE = 10; DistributedVectorFactory factory(PROBLEM_SIZE); // Source vector = [-1 1 -2 2 ... -PROBLEM_SIZE/2+1 PROBLEM_SIZE/2+1] Vec interleaved_vec=factory.CreateVec(2); DistributedVector interleaved_dist_vec = factory.CreateDistributedVector(interleaved_vec); DistributedVector::Stripe first_variable(interleaved_dist_vec, 0); DistributedVector::Stripe second_variable(interleaved_dist_vec, 1); for (DistributedVector::Iterator index = interleaved_dist_vec.Begin(); index!= interleaved_dist_vec.End(); ++index) { first_variable[index] = -1.0 * (index.Global+1); second_variable[index] = (index.Global+1); } // Destination vectors. It is important they have a compatible parallel layout with the source vector, hence the 2nd param of CreateVec() PetscInt local_num_rows; VecGetLocalSize(interleaved_vec, &local_num_rows); Vec first_variable_vec = PetscTools::CreateVec(PROBLEM_SIZE, local_num_rows/2); Vec second_variable_vec = PetscTools::CreateVec(PROBLEM_SIZE, local_num_rows/2); // Setup and perform scatter operation VecScatter first_variable_context; VecScatter second_variable_context; PetscVecTools::SetupInterleavedVectorScatterGather(interleaved_vec, first_variable_context, second_variable_context); #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 6) //PETSc 3.6 or later // When used in the preconditioner code within a KSP solve, the main (bidomain) vector will be locked. // Lock the vector so that modern PETSc (3.6) won't want to change it VecLockPush(interleaved_vec); #endif PetscVecTools::DoInterleavedVecScatter(interleaved_vec, first_variable_context, first_variable_vec, second_variable_context, second_variable_vec); #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 6) //PETSc 3.6 or later // Unlock the vector (for symmetry) VecLockPop(interleaved_vec); #endif // Check destination vectors are [-1 -2 -3 ...] and [1 2 3 ...] respectively. DistributedVector dist_1st_var_vec = factory.CreateDistributedVector(first_variable_vec); DistributedVector dist_2nd_var_vec = factory.CreateDistributedVector(second_variable_vec); for (DistributedVector::Iterator index = dist_1st_var_vec.Begin(); index!= dist_1st_var_vec.End(); ++index) { TS_ASSERT_EQUALS(dist_1st_var_vec[index], -1.0 * (index.Global+1)); TS_ASSERT_EQUALS(dist_2nd_var_vec[index], index.Global+1); } PetscTools::Destroy(interleaved_vec); PetscTools::Destroy(first_variable_vec); PetscTools::Destroy(second_variable_vec); VecScatterDestroy(PETSC_DESTROY_PARAM(first_variable_context)); VecScatterDestroy(PETSC_DESTROY_PARAM(second_variable_context)); }
void TestInterleavedVecGather() { // Vectors will be twice PROBLEM_SIZE, since this is to be used in bidomain code. const unsigned PROBLEM_SIZE = 10; DistributedVectorFactory factory(PROBLEM_SIZE); // Destination vector Vec interleaved_vec=factory.CreateVec(2); // Source vectors. It is important they have a compatible parallel layout with the destination vector, hence the 2nd param of CreateVec() PetscInt local_num_rows; VecGetLocalSize(interleaved_vec, &local_num_rows); Vec first_variable_vec = PetscTools::CreateVec(PROBLEM_SIZE, local_num_rows/2); Vec second_variable_vec = PetscTools::CreateVec(PROBLEM_SIZE, local_num_rows/2); // Fill in source vectors. DistributedVector dist_first_variable_vec = factory.CreateDistributedVector(first_variable_vec); DistributedVector dist_second_variable_vec = factory.CreateDistributedVector(second_variable_vec); for (DistributedVector::Iterator index = dist_first_variable_vec.Begin(); index!= dist_first_variable_vec.End(); ++index) { dist_first_variable_vec[index] = -1.0 * (index.Global+1); dist_second_variable_vec[index] = (index.Global+1); } // Setup and perform gather operation VecScatter first_variable_context; VecScatter second_variable_context; PetscVecTools::SetupInterleavedVectorScatterGather(interleaved_vec, first_variable_context, second_variable_context); PetscVecTools::DoInterleavedVecGather(interleaved_vec, first_variable_context, first_variable_vec, second_variable_context, second_variable_vec); // Check that the destination vector looks like [ -1 1 -2 2 ... -PROBLEM_SIZE/2+1 PROBLEM_SIZE/2+1 ] DistributedVector dist_interleaved_vec = factory.CreateDistributedVector(interleaved_vec); DistributedVector::Stripe dist_inter_vec_1st_var(dist_interleaved_vec, 0); DistributedVector::Stripe dist_inter_vec_2nd_var(dist_interleaved_vec, 1); for (DistributedVector::Iterator index = dist_interleaved_vec.Begin(); index!= dist_interleaved_vec.End(); ++index) { TS_ASSERT_EQUALS(dist_inter_vec_1st_var[index], -1.0 * (index.Global+1)); TS_ASSERT_EQUALS(dist_inter_vec_2nd_var[index], index.Global+1 ); } PetscTools::Destroy(interleaved_vec); PetscTools::Destroy(first_variable_vec); PetscTools::Destroy(second_variable_vec); VecScatterDestroy(PETSC_DESTROY_PARAM(first_variable_context)); VecScatterDestroy(PETSC_DESTROY_PARAM(second_variable_context)); }
VentilationProblem::~VentilationProblem() { /* Remove the PETSc context used in the iterative solver */ if (mTerminalInteractionMatrix) { PetscTools::Destroy(mTerminalInteractionMatrix); PetscTools::Destroy(mTerminalFluxChangeVector); PetscTools::Destroy(mTerminalPressureChangeVector); KSPDestroy(PETSC_DESTROY_PARAM(mTerminalKspSolver)); } }
void ReplicatableVector::RemovePetscContext() { if (mToAll != NULL) { VecScatterDestroy(PETSC_DESTROY_PARAM(mToAll)); mToAll = NULL; } if (mReplicated != NULL) { PetscTools::Destroy(mReplicated); mReplicated = NULL; } if (mpData != NULL) { delete[] mpData; mpData = NULL; } }
void VentilationProblem::SolveFromPressureWithSnes() { assert( !mFluxGivenAtInflow ); // It's not a direct solve if (mTerminalInteractionMatrix == NULL) { SetupIterativeSolver(); } SNES snes; SNESCreate(PETSC_COMM_SELF, &snes); // Set the residual creation function (direct solve flux->pressure followed by pressure matching) SNESSetFunction(snes, mTerminalPressureChangeVector /*residual*/ , &ComputeSnesResidual, this); // The approximate Jacobian has been precomputed so we are going to wing it SNESSetJacobian(snes, mTerminalInteractionMatrix, mTerminalInteractionMatrix, /*&ComputeSnesJacobian*/ NULL, this); #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x SNESSetLagJacobian(snes, -1 /*Never rebuild Jacobian*/); #else SNESSetType(snes, SNESLS); #endif #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later SNESSetType(snes, SNESNEWTONLS); #else SNESSetType(snes, SNESLS); #endif // Set the relative tolerance on the residual // Also set the absolute tolerance - useful for when the solver is started from the correct answer SNESSetTolerances(snes, 1.0e-16/*abs_tol*/, 1e-15/*r_tol*/, PETSC_DEFAULT/*step_tol*/, PETSC_DEFAULT, PETSC_DEFAULT); #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 3) //PETSc 3.3 SNESLineSearch linesearch; SNESGetSNESLineSearch(snes, &linesearch); SNESLineSearchSetType(linesearch, "bt"); //Use backtracking search as default #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later SNESLineSearch linesearch; SNESGetLineSearch(snes, &linesearch); SNESLineSearchSetType(linesearch, "bt"); //Use backtracking search as default #endif SNESSetFromOptions(snes); #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5) // Seems to want the preconditioner to be explicitly set to none now // Copied this from the similar PETSc example at: // http://www.mcs.anl.gov/petsc/petsc-current/src/snes/examples/tutorials/ex1.c // Which got it to work... KSP ksp; SNESGetKSP(snes,&ksp); PC pc; KSPGetPC(ksp,&pc); PCSetType(pc,PCNONE); #endif #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 SNESSolve(snes, mTerminalFluxChangeVector); #else SNESSolve(snes, PETSC_NULL, mTerminalFluxChangeVector); #endif ///\todo #2300 If used with time-stepping we should maintain a permanent SNES object SNESDestroy(PETSC_DESTROY_PARAM(snes)); }
void PetscMatTools::ZeroRowsWithValueOnDiagonal(Mat matrix, std::vector<unsigned>& rRows, double diagonalValue) { Finalise(matrix); /* * Important! PETSc by default will destroy the sparsity structure for this row and deallocate memory * when the row is zeroed, and if there is a next timestep, the memory will have to reallocated when * assembly to done again. This can kill performance. The following makes sure the zeroed rows are kept. * * Note: if the following lines are called, then once MatZeroRows() is called below, there will be an * error if some rows have no entries at all. Hence for problems with dummy variables, Stokes flow for * example, the identity block needs to be added before dirichlet BCs. */ #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later MatSetOption(matrix, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE); #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 0) //PETSc 3.0 MatSetOption(matrix, MAT_KEEP_ZEROED_ROWS, PETSC_TRUE); #else //PETSc 2.x.x MatSetOption(matrix, MAT_KEEP_ZEROED_ROWS); #endif PetscInt* rows = new PetscInt[rRows.size()]; for (unsigned i=0; i<rRows.size(); i++) { rows[i] = rRows[i]; } #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2 IS is; ISCreateGeneral(PETSC_COMM_WORLD, rRows.size(), rows, &is); MatZeroRows(matrix, is, &diagonalValue); ISDestroy(PETSC_DESTROY_PARAM(is)); /* * [2]PETSC ERROR: MatMissingDiagonal_SeqAIJ() line 1011 in src/mat/impls/aij/seq/aij.c [2]PETSC ERROR: PETSc has generated inconsistent data! [2]PETSC ERROR: Matrix is missing diagonal number 15! [2]PETSC ERROR: MatILUFactorSymbolic_SeqAIJ() line 906 in src/mat/impls/aij/seq/aijfact.c * * // There appears to be a problem with MatZeroRows not setting diagonals correctly // While we are supporting PETSc 2.2, we have to do this the slow way AssembleFinal(matrix); PetscInt lo, hi; GetOwnershipRange(matrix, lo, hi); PetscInt size=GetSize(matrix); ///assert(rRows.size() == 1); for (unsigned index=0; index<rRows.size(); index++) { PetscInt row = rRows[index]; if (row >= lo && row < hi) { std::vector<unsigned> non_zero_cols; // This row is local, so zero it. for (PetscInt column = 0; column < size; column++) { if (GetElement(matrix, row, column) != 0.0) { non_zero_cols.push_back(column); } } // Separate "gets" from "sets" so that we don't have to keep going into "assembled" mode for (unsigned i=0; i<non_zero_cols.size();i++) { SetElement(matrix, row, non_zero_cols[i], 0.0); } // Set the diagonal SetElement(matrix, row, row, diagonalValue); } // Everyone communicate after row is finished AssembleFinal(matrix); } */ #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 2) //PETSc 3.2 or later MatZeroRows(matrix, rRows.size(), rows, diagonalValue , NULL, NULL); #else MatZeroRows(matrix, rRows.size(), rows, diagonalValue); #endif delete [] rows; }