Class ae108::solve::GeneralizedTAOSolver
template <class Assembler>
ClassList > ae108 > solve > GeneralizedTAOSolver
Classes
Type | Name |
---|---|
struct | BoundaryConditionContainer |
Public Types
Type | Name |
---|---|
typedef std::function< void(const cpppetsc::distributed< vector_type > &, double, real_type *)> | DistributedEnergyAssembler |
typedef std::function< void(const cpppetsc::distributed< vector_type > &, double, cpppetsc::distributed< vector_type > *)> | DistributedForceVectorAssembler |
typedef std::function< void(const cpppetsc::distributed< vector_type > &, double, matrix_type *)> | DistributedStiffnessMatrixAssembler |
typedef std::function< void(const cpppetsc::local< vector_type > &, double, real_type *)> | LocalEnergyAssembler |
typedef std::function< void(const cpppetsc::local< vector_type > &, double, cpppetsc::local< vector_type > *)> | LocalForceVectorAssembler |
typedef std::function< void(const cpppetsc::local< vector_type > &, double, matrix_type *)> | LocalStiffnessMatrixAssembler |
typedef typename mesh_type::matrix_type | matrix_type |
typedef typename assembly::MeshTypeTrait< Assembler >::type | mesh_type |
typedef typename assembly::PolicyTypeTrait< Assembler >::type | policy_type |
typedef typename mesh_type::real_type | real_type |
typedef typename mesh_type::size_type | size_type |
typedef typename mesh_type::value_type | value_type |
typedef typename mesh_type::vector_type | vector_type |
Public Functions
Type | Name |
---|---|
GeneralizedTAOSolver (const mesh_type * mesh) |
|
cpppetsc::distributed< vector_type > | computeSolution (const BoundaryConditionContainer & boundaryConditions, cpppetsc::distributed< vector_type > initialGuess, const double time, const Assembler *const assembler) const Calls computeSolution with the default assembler calls (ie. calls assembleEnergy, assembleForceVector or assembleStiffnessMatrix). Passes on the rest of the arguments. |
cpppetsc::distributed< vector_type > | computeSolution (const BoundaryConditionContainer & boundaryConditions, cpppetsc::distributed< vector_type > initialGuess, const double time, LocalEnergyAssembler assembleEnergy, LocalForceVectorAssembler assembleForceVector, LocalStiffnessMatrixAssembler assembleStiffnessMatrix) const Calls computeSolution wrapping up the local assembler calls into distributed assembler calls. Passes on the rest of the arguments. |
cpppetsc::distributed< vector_type > | computeSolution (const BoundaryConditionContainer & boundaryConditions, cpppetsc::distributed< vector_type > initialGuess, const double time, DistributedEnergyAssembler assembleEnergy, DistributedForceVectorAssembler assembleForceVector, DistributedStiffnessMatrixAssembler assembleStiffnessMatrix) const The goal of this function is to solve E(x) = min under the condition that the residual is equal to zero. |
Public Types Documentation
typedef DistributedEnergyAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::DistributedEnergyAssembler = std::function<void( const cpppetsc::distributed<vector_type> &, double, real_type *)>;
typedef DistributedForceVectorAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::DistributedForceVectorAssembler = std::function<void(const cpppetsc::distributed<vector_type> &, double, cpppetsc::distributed<vector_type> *)>;
typedef DistributedStiffnessMatrixAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::DistributedStiffnessMatrixAssembler = std::function<void( const cpppetsc::distributed<vector_type> &, double, matrix_type *)>;
typedef LocalEnergyAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::LocalEnergyAssembler = std::function<void( const cpppetsc::local<vector_type> &, double, real_type *)>;
typedef LocalForceVectorAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::LocalForceVectorAssembler = std::function<void(const cpppetsc::local<vector_type> &, double, cpppetsc::local<vector_type> *)>;
typedef LocalStiffnessMatrixAssembler
using ae108::solve::GeneralizedTAOSolver< Assembler >::LocalStiffnessMatrixAssembler = std::function<void( const cpppetsc::local<vector_type> &, double, matrix_type *)>;
typedef matrix_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::matrix_type = typename mesh_type::matrix_type;
typedef mesh_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::mesh_type = typename assembly::MeshTypeTrait<Assembler>::type;
typedef policy_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::policy_type = typename assembly::PolicyTypeTrait<Assembler>::type;
typedef real_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::real_type = typename mesh_type::real_type;
typedef size_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::size_type = typename mesh_type::size_type;
typedef value_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::value_type = typename mesh_type::value_type;
typedef vector_type
using ae108::solve::GeneralizedTAOSolver< Assembler >::vector_type = typename mesh_type::vector_type;
Public Functions Documentation
function GeneralizedTAOSolver
explicit ae108::solve::GeneralizedTAOSolver::GeneralizedTAOSolver (
const mesh_type * mesh
)
Parameters:
mesh
A valid nonzero pointer to a mesh_type instance.
function computeSolution [1/3]
Calls computeSolution with the default assembler calls (ie. calls assembleEnergy, assembleForceVector or assembleStiffnessMatrix). Passes on the rest of the arguments.
cpppetsc::distributed < vector_type > ae108::solve::GeneralizedTAOSolver::computeSolution (
const BoundaryConditionContainer & boundaryConditions,
cpppetsc::distributed < vector_type > initialGuess,
const double time,
const Assembler *const assembler
) const
Parameters:
boundaryConditions
The nonlinear boundary conditions to apply.initialGuess
The global vector to start iterating at.time
This will be used to configure the assembler.assembler
Valid nonzero pointer.
function computeSolution [2/3]
Calls computeSolution wrapping up the local assembler calls into distributed assembler calls. Passes on the rest of the arguments.
cpppetsc::distributed < vector_type > ae108::solve::GeneralizedTAOSolver::computeSolution (
const BoundaryConditionContainer & boundaryConditions,
cpppetsc::distributed < vector_type > initialGuess,
const double time,
LocalEnergyAssembler assembleEnergy,
LocalForceVectorAssembler assembleForceVector,
LocalStiffnessMatrixAssembler assembleStiffnessMatrix
) const
Parameters:
boundaryConditions
The nonlinear boundary conditions to apply. Note that empty functions (residual, jacobian) are interpreted as functions that do not change the output.initialGuess
The global vector to start iterating at.time
This will be used to configure the assembler.assembleEnergy
A valid callable. It will be called to assemble the local energy.assembleForceVector
A valid callable. It will be called to assemble the local force vector.assembleStiffnessMatrix
A valid callable. It will be called to assemble the stiffness matrix.
function computeSolution [3/3]
The goal of this function is to solve E(x) = min under the condition that the residual is equal to zero.
cpppetsc::distributed < vector_type > ae108::solve::GeneralizedTAOSolver::computeSolution (
const BoundaryConditionContainer & boundaryConditions,
cpppetsc::distributed < vector_type > initialGuess,
const double time,
DistributedEnergyAssembler assembleEnergy,
DistributedForceVectorAssembler assembleForceVector,
DistributedStiffnessMatrixAssembler assembleStiffnessMatrix
) const
Remark:
Does not update internal variables.
Parameters:
boundaryConditions
The nonlinear boundary conditions to apply.initialGuess
The global vector to start iterating at.time
This will be used to configure the assembler.assembleEnergy
A valid callable. It will be called to assemble the distributed energy.assembleForceVector
A valid callable. It will be called to assemble the distributed force vector.assembleStiffnessMatrix
A valid callable. It will be called to assemble the stiffness matrix.
The documentation for this class was generated from the following file solve/src/include/ae108/solve/GeneralizedTAOSolver.h