99 INTEGER(INTG),
INTENT(IN) :: SOLUTION_METHOD
100 INTEGER(INTG),
INTENT(OUT) :: Err
105 enters(
"FSI_EQUATIONS_SET_SOLUTION_METHOD_SET",err,error,*999)
107 IF(
ASSOCIATED(equations_set))
THEN 108 IF(.NOT.
ALLOCATED(equations_set%SPECIFICATION))
THEN 109 CALL flagerror(
"Equations set specification is not allocated.",err,error,*999)
110 ELSE IF(
SIZE(equations_set%SPECIFICATION,1)/=3)
THEN 111 CALL flagerror(
"Equations set specification must have three entries for a "// &
112 &
"finite elasticity Navier-Stokes class equations set.",err,error,*999)
114 SELECT CASE(equations_set%SPECIFICATION(3))
116 SELECT CASE(solution_method)
120 CALL flagerror(
"Not implemented.",err,error,*999)
122 CALL flagerror(
"Not implemented.",err,error,*999)
124 CALL flagerror(
"Not implemented.",err,error,*999)
126 CALL flagerror(
"Not implemented.",err,error,*999)
128 CALL flagerror(
"Not implemented.",err,error,*999)
130 local_error=
"The specified solution method of "//
trim(
number_to_vstring(solution_method,
"*",err,error))//
" is invalid." 131 CALL flagerror(local_error,err,error,*999)
134 local_error=
"Equations set subtype of "//
trim(
number_to_vstring(equations_set%SPECIFICATION(3),
"*",err,error))// &
135 &
" is not valid for a finite elasticity navier stokes equation type of a multi physics equations set class." 136 CALL flagerror(local_error,err,error,*999)
139 CALL flagerror(
"Equations set is not associated.",err,error,*999)
142 exits(
"FSI_EQUATIONS_SET_SOLUTION_METHOD_SET")
144 999 errorsexits(
"FSI_EQUATIONS_SET_SOLUTION_METHOD_SET",err,error)
158 INTEGER(INTG),
INTENT(OUT) :: Err
162 enters(
"FSI_EQUATIONS_SET_SETUP",err,error,*999)
164 CALL flagerror(
"FSI_EQUATIONS_SET_SETUP is not implemented.",err,error,*999)
166 exits(
"FSI_EQUATIONS_SET_SETUP")
168 999 errorsexits(
"FSI_EQUATIONS_SET_SETUP",err,error)
181 INTEGER(INTG),
INTENT(IN) :: ELEMENT_NUMBER
182 INTEGER(INTG),
INTENT(OUT) :: Err
185 enters(
"FSI_FINITE_ELEMENT_CALCULATE",err,error,*999)
187 CALL flagerror(
"FSI_FINITE_ELEMENT_CALCULATE is not implemented.",err,error,*999)
189 exits(
"FSI_FINITE_ELEMENT_CALCULATE")
191 999 errorsexits(
"FSI_FINITE_ELEMENT_CALCULATE",err,error)
204 INTEGER(INTG),
INTENT(IN) :: specification(:)
205 INTEGER(INTG),
INTENT(OUT) :: err
209 enters(
"FSI_EquationsSetSpecificationSet",err,error,*999)
211 CALL flagerror(
"FSI_EquationsSetSpecificationSet is not implemented.",err,error,*999)
213 exits(
"FSI_EquationsSetSpecificationSet")
215 999
errors(
"FSI_EquationsSetSpecificationSet",err,error)
216 exits(
"FSI_EquationsSetSpecificationSet")
230 INTEGER(INTG),
INTENT(IN) :: problemSpecification(:)
231 INTEGER(INTG),
INTENT(OUT) :: err
235 INTEGER(INTG) :: problemSubtype
237 enters(
"FSI_ProblemSpecificationSet",err,error,*999)
239 IF(
ASSOCIATED(problem))
THEN 240 IF(
SIZE(problemspecification,1)==3)
THEN 241 problemsubtype=problemspecification(3)
242 SELECT CASE(problemsubtype)
246 localerror=
"The third problem specification of "//
trim(
numbertovstring(problemsubtype,
"*",err,error))// &
247 &
" is not valid for a finite elasticity Navier-Stokes type of a multi physics problem." 248 CALL flagerror(localerror,err,error,*999)
250 IF(
ALLOCATED(problem%specification))
THEN 251 CALL flagerror(
"Problem specification is already allocated.",err,error,*999)
253 ALLOCATE(problem%specification(3),stat=err)
254 IF(err/=0)
CALL flagerror(
"Could not allocate problem specification.",err,error,*999)
260 CALL flagerror(
"Finite elasticity Navier-Stokes problem specificaion must have three entries.",err,error,*999)
263 CALL flagerror(
"Problem is not associated.",err,error,*999)
266 exits(
"FSI_ProblemSpecificationSet")
268 999
errors(
"FSI_ProblemSpecificationSet",err,error)
269 exits(
"FSI_ProblemSpecificationSet")
284 INTEGER(INTG),
INTENT(OUT) :: Err
289 TYPE(
solver_type),
POINTER :: SOLVER,MOVING_MESH_SOLVER
290 TYPE(
solver_equations_type),
POINTER :: SOLVER_EQUATIONS,MovingMeshSolverEquations,MOVING_MESH_SOLVER_EQUATIONS
291 TYPE(
solvers_type),
POINTER :: MovingMeshSolvers,SOLVERS
294 enters(
"FSI_PROBLEM_SETUP",err,error,*999)
296 NULLIFY(control_loop)
297 NULLIFY(movingmeshsubloop)
298 NULLIFY(elasticitysubloop)
300 NULLIFY(moving_mesh_solver)
301 NULLIFY(solver_equations)
302 NULLIFY(moving_mesh_solver_equations)
303 NULLIFY(movingmeshsolverequations)
305 NULLIFY(movingmeshsolvers)
307 IF(
ASSOCIATED(problem))
THEN 308 IF(
ALLOCATED(problem%specification))
THEN 309 IF(.NOT.
ALLOCATED(problem%specification))
THEN 310 CALL flagerror(
"Problem specification is not allocated.",err,error,*999)
311 ELSE IF(
SIZE(problem%specification,1)<3)
THEN 312 CALL flagerror(
"Problem specification must have three entries for a finite elasticity-Darcy problem.", &
316 CALL flagerror(
"Problem specification is not allocated.",err,error,*999)
318 SELECT CASE(problem%SPECIFICATION(3))
321 SELECT CASE(problem_setup%SETUP_TYPE)
323 SELECT CASE(problem_setup%ACTION_TYPE)
329 local_error=
"The action type of "//
trim(
number_to_vstring(problem_setup%ACTION_TYPE,
"*",err,error))// &
331 &
" is invalid for an finite elasticity ALE navier stokes equation." 332 CALL flagerror(local_error,err,error,*999)
335 SELECT CASE(problem_setup%ACTION_TYPE)
343 control_loop_root=>problem%CONTROL_LOOP
347 local_error=
"The action type of "//
trim(
number_to_vstring(problem_setup%ACTION_TYPE,
"*",err,error))// &
349 &
" is invalid for a finite elasticity navier stokes equation." 350 CALL flagerror(local_error,err,error,*999)
354 control_loop_root=>problem%CONTROL_LOOP
356 SELECT CASE(problem_setup%ACTION_TYPE)
381 local_error=
"The action type of "//
trim(
number_to_vstring(problem_setup%ACTION_TYPE,
"*",err,error))// &
383 &
" is invalid for a finite elasticity navier stokes equation." 384 CALL flagerror(local_error,err,error,*999)
387 SELECT CASE(problem_setup%ACTION_TYPE)
390 control_loop_root=>problem%CONTROL_LOOP
410 control_loop_root=>problem%CONTROL_LOOP
424 local_error=
"The action type of "//
trim(
number_to_vstring(problem_setup%ACTION_TYPE,
"*",err,error))// &
426 &
" is invalid for a finite elasticity navier stokes equation." 427 CALL flagerror(local_error,err,error,*999)
430 local_error=
"The setup type of "//
trim(
number_to_vstring(problem_setup%SETUP_TYPE,
"*",err,error))// &
431 &
" is invalid for a finite elasticity ALE navier stokes equation." 432 CALL flagerror(local_error,err,error,*999)
435 local_error=
"The problem subtype of "//
trim(
number_to_vstring(problem%SPECIFICATION(3),
"*",err,error))// &
436 &
" does not equal a standard finite elasticity navier stokes equation subtype." 437 CALL flagerror(local_error,err,error,*999)
440 CALL flagerror(
"Problem is not associated.",err,error,*999)
443 exits(
"FSI_PROBLEM_SETUP")
445 999 errorsexits(
"FSI_PROBLEM_SETUP",err,error)
459 INTEGER(INTG),
INTENT(OUT) :: Err
465 enters(
"FSI_PRE_SOLVE",err,error,*999)
467 IF(
ASSOCIATED(controlloop))
THEN 468 IF(
ASSOCIATED(solver))
THEN 469 IF(
ASSOCIATED(controlloop%problem))
THEN 470 IF(.NOT.
ALLOCATED(controlloop%problem%specification))
THEN 471 CALL flagerror(
"Problem specification is not allocated.",err,error,*999)
472 ELSE IF(
SIZE(controlloop%problem%specification,1)<3)
THEN 473 CALL flagerror(
"Problem specification must have three entries for an elasticity Navier-Stokes problem.",err,error,*999)
475 SELECT CASE(controlloop%problem%specification(3))
480 CALL navier_stokes_pre_solve(solver,err,error,*999)
484 CALL flagerror(
"Incorrect loop type. Must be time loop.",err,error,*999)
487 local_error=
"Problem subtype "//
trim(
number_to_vstring(controlloop%PROBLEM%SPECIFICATION(3),
"*",err,error))// &
488 &
" is not valid for a finite elasticity navier stokes type of a multi physics problem class." 489 CALL flagerror(local_error,err,error,*999)
492 CALL flagerror(
"Problem is not associated.",err,error,*999)
495 CALL flagerror(
"Solver is not associated.",err,error,*999)
498 CALL flagerror(
"Control loop is not associated.",err,error,*999)
501 exits(
"FSI_PRE_SOLVE")
503 999 errorsexits(
"FSI_PRE_SOLVE",err,error)
517 INTEGER(INTG),
INTENT(OUT) :: Err
523 enters(
"FSI_POST_SOLVE",err,error,*999)
527 IF(
ASSOCIATED(controlloop))
THEN 528 IF(
ASSOCIATED(solver))
THEN 529 IF(
ASSOCIATED(controlloop%problem))
THEN 530 IF(.NOT.
ALLOCATED(controlloop%problem%specification))
THEN 531 CALL flagerror(
"Problem specification is not allocated.",err,error,*999)
532 ELSE IF(
SIZE(controlloop%problem%specification,1)<3)
THEN 533 CALL flagerror(
"Problem specification must have three entries for an elasticity Navier-Stokes problem.",err,error,*999)
535 SELECT CASE(controlloop%PROBLEM%SPECIFICATION(3))
542 IF(
ASSOCIATED(solver2%DYNAMIC_SOLVER))
THEN 543 solver2%DYNAMIC_SOLVER%ALE=.true.
545 CALL flagerror(
"Dynamic solver is not associated for ALE problem.",err,error,*999)
570 local_error=
"Problem subtype "//
trim(
number_to_vstring(controlloop%PROBLEM%SPECIFICATION(3),
"*",err,error))// &
571 &
" for a FiniteElasticity-NavierStokes type of a multi physics problem class has unknown solver solve type." 572 CALL flagerror(local_error,err,error,*999)
575 local_error=
"Problem subtype "//
trim(
number_to_vstring(controlloop%PROBLEM%SPECIFICATION(3),
"*",err,error))// &
576 &
" is not valid for a FiniteElasticity-NavierStokes type of a multi physics problem class." 577 CALL flagerror(local_error,err,error,*999)
580 CALL flagerror(
"Problem is not associated.",err,error,*999)
583 CALL flagerror(
"Solver is not associated.",err,error,*999)
586 CALL flagerror(
"Control loop is not associated.",err,error,*999)
589 exits(
"FSI_POST_SOLVE")
591 999 errorsexits(
"FSI_POST_SOLVE",err,error)
604 INTEGER(INTG),
INTENT(OUT) :: Err
609 enters(
"FSI_CONTROL_LOOP_PRE_LOOP",err,error,*999)
611 CALL flagerror(
"FSI_CONTROL_LOOP_PRE_LOOP not implemented.",err,error,*999)
613 exits(
"FSI_CONTROL_LOOP_PRE_LOOP")
615 999 errorsexits(
"FSI_CONTROL_LOOP_PRE_LOOP",err,error)
628 INTEGER(INTG),
INTENT(OUT) :: Err
632 TYPE(
solver_type),
POINTER :: DynamicSolver,LinearSolver
637 TYPE(
field_type),
POINTER :: SolidGeometricField,InterfaceGeometricField,SolidDependentField
641 TYPE(
varying_string) :: Method,SolidFileName,FluidFileName,InterfaceFileName
642 REAL(DP) :: StartTime,CurrentTime,TimeIncrement,TimeStepNumber,Value
643 INTEGER(INTG) :: EquationsSetIndex,InterfaceNodeNumber,InterfaceNodeComponent,NumberOfComponents
644 LOGICAL :: FluidEquationsSetFound,SolidEquationsSetFound=.false.
646 enters(
"FSI_CONTROL_LOOP_POST_LOOP",err,error,*999)
648 NULLIFY(dynamicsolver)
649 NULLIFY(linearsolver)
651 NULLIFY(dynamicsolverequations)
652 NULLIFY(solidequationsset)
653 NULLIFY(fluidequationsset)
654 NULLIFY(solidgeometricfield)
655 NULLIFY(interfacegeometricfield)
656 NULLIFY(soliddependentfield)
657 NULLIFY(interfacecondition)
659 NULLIFY(interfacenodes)
662 IF(.NOT.
ASSOCIATED(controlloop))
CALL flagerror(
"Main control loop not associated.",err,error,*999)
663 IF(.NOT.
ASSOCIATED(controlloop%SOLVERS))
CALL flagerror(
"Solvers are not associated.",err,error,*999)
667 timeloop=>controlloop%TIME_LOOP
668 IF(.NOT.
ASSOCIATED(timeloop))
CALL flagerror(
"Time loop not associated.",err,error,*999)
670 starttime=timeloop%START_TIME
671 currenttime=timeloop%CURRENT_TIME
672 timeincrement=timeloop%TIME_INCREMENT
673 timestepnumber=(currenttime-starttime)/timeincrement
676 CALL navier_stokes_pre_solve_ale_update_mesh(dynamicsolver,err,error,*999)
679 dynamicsolverequations=>dynamicsolver%SOLVER_EQUATIONS
680 IF(
ASSOCIATED(dynamicsolverequations))
THEN 681 dynamicsolvermapping=>dynamicsolverequations%SOLVER_MAPPING
682 IF(
ASSOCIATED(dynamicsolvermapping))
THEN 684 fluidequationssetfound=.false.
685 solidequationssetfound=.false.
686 DO WHILE((equationssetindex<=dynamicsolvermapping%NUMBER_OF_EQUATIONS_SETS &
687 & .AND..NOT.solidequationssetfound) &
688 & .OR.(equationssetindex<=dynamicsolvermapping%NUMBER_OF_EQUATIONS_SETS &
689 & .AND..NOT.fluidequationssetfound))
690 equationsset=>dynamicsolvermapping%EQUATIONS_SETS(equationssetindex)%PTR
695 solidequationsset=>equationsset
696 solidequationssetfound=.true.
700 fluidequationsset=>equationsset
701 fluidequationssetfound=.true.
703 CALL flagerror(
"Invalid equations sets associated with dynamic solver for FSI.", err,error,*999)
705 equationssetindex=equationssetindex+1
707 IF(.NOT.solidequationssetfound)
CALL flagerror(
"Could not find solid equations set for FSI.",err,error,*999)
708 IF(.NOT.fluidequationssetfound)
CALL flagerror(
"Could not find fluid equations set for FSI.",err,error,*999)
709 solidgeometricfield=>solidequationsset%GEOMETRY%GEOMETRIC_FIELD
710 IF(
ASSOCIATED(solidgeometricfield))
THEN 711 CALL field_number_of_components_get(solidgeometricfield,field_u_variable_type,numberofcomponents,err,error,*999)
712 IF(dynamicsolvermapping%NUMBER_OF_INTERFACE_CONDITIONS>1)
THEN 713 CALL flagerror(
"Invalid number of interface conditions. Must be 1 for FSI.",err,error,*999)
715 soliddependentfield=>solidequationsset%DEPENDENT%DEPENDENT_FIELD
716 IF(
ASSOCIATED(soliddependentfield))
THEN 717 interfacecondition=>dynamicsolvermapping%INTERFACE_CONDITIONS(1)%PTR
718 IF(
ASSOCIATED(interfacecondition))
THEN 719 fsinterface=>interfacecondition%INTERFACE
720 IF(
ASSOCIATED(fsinterface))
THEN 721 interfacenodes=>fsinterface%NODES
722 IF(
ASSOCIATED(interfacenodes))
THEN 723 interfacegeometricfield=>interfacecondition%GEOMETRY%GEOMETRIC_FIELD
724 IF(
ASSOCIATED(interfacegeometricfield))
THEN 727 DO interfacenodenumber=1,interfacenodes%NUMBER_OF_NODES
728 DO interfacenodecomponent=1,numberofcomponents
730 CALL field_parameter_set_get_node(soliddependentfield,field_u_variable_type,field_values_set_type, &
731 & 1,1,interfacenodes%COUPLED_NODES(1,interfacenodenumber),interfacenodecomponent,
Value, &
733 CALL field_parameter_set_update_node(interfacegeometricfield,field_u_variable_type, &
734 & field_values_set_type,1,1,interfacenodenumber,interfacenodecomponent,
Value,err,error,*999)
737 CALL field_parameter_set_update_start(interfacegeometricfield, &
738 & field_u_variable_type,field_values_set_type,err,error,*999)
739 CALL field_parameter_set_update_finish(interfacegeometricfield, &
740 & field_u_variable_type,field_values_set_type,err,error,*999)
745 interfacefilename=
"./output/Interface/Interface"//
trim(
number_to_vstring(int(timestepnumber),
"*",err,error))
748 IF(.NOT.
ASSOCIATED(solidequationsset%REGION))
CALL flagerror(
"Solid region not associated.", &
750 IF(.NOT.
ASSOCIATED(solidequationsset%REGION%FIELDS))
CALL flagerror(
"Solid fields not associated.", &
757 IF(.NOT.
ASSOCIATED(fluidequationsset%REGION))
CALL flagerror(
"Fluid region not associated.", &
759 IF(.NOT.
ASSOCIATED(fluidequationsset%REGION%FIELDS))
CALL flagerror(
"Fluid fields not associated.", &
765 IF(.NOT.
ASSOCIATED(fsinterface%FIELDS))
CALL flagerror(
"Interface fields not associated.",err,error,*999)
772 CALL flagerror(
"Interface geometric field not associated.",err,error,*999)
775 CALL flagerror(
"Interface nodes not associated.",err,error,*999)
778 CALL flagerror(
"Interface not associated.",err,error,*999)
781 CALL flagerror(
"Interface condition not associated.",err,error,*999)
784 CALL flagerror(
"Solid dependent field not associated.",err,error,*999)
787 CALL flagerror(
"Solid geometric field not associated.",err,error,*999)
790 CALL flagerror(
"Dynamic solver mapping not associated.",err,error,*999)
793 CALL flagerror(
"Dynamic solver equations not associated.",err,error,*999)
796 exits(
"FSI_CONTROL_LOOP_POST_LOOP")
798 999 errorsexits(
"FSI_CONTROL_LOOP_POST_LOOP",err,error)
integer(intg), parameter equations_set_fem_solution_method
Finite Element Method solution method.
This module contains all basis function routines.
subroutine, public enters(NAME, ERR, ERROR,)
Records the entry into the named procedure and initialises the error code.
subroutine, public solvers_create_finish(SOLVERS, ERR, ERROR,)
Finish the creation of solvers.
integer(intg), parameter, public control_loop_progress_output
Progress output from control loop.
integer(intg), parameter equations_set_gfem_solution_method
Grid-based Finite Element Method solution method.
integer(intg), parameter problem_control_time_loop_type
Time control loop.
Contains information on a time iteration control loop.
integer(intg), parameter problem_setup_control_type
Solver setup for a problem.
This module handles all problem wide constants.
integer(intg), parameter solver_equations_first_order_dynamic
Solver equations are first order dynamic.
integer(intg), parameter, public control_loop_node
The identifier for a each "leaf" node in a control loop.
subroutine, public solver_dynamic_order_set(SOLVER, ORDER, ERR, ERROR,)
Sets/changes the order for a dynamic solver.
Converts a number to its equivalent varying string representation.
Contains information on the type of solver to be used.
integer(intg), parameter, public solver_petsc_library
PETSc solver library.
subroutine, public solvers_number_set(SOLVERS, NUMBER_OF_SOLVERS, ERR, ERROR,)
Sets/changes the number of solvers.
integer(intg), parameter, public solver_dynamic_crank_nicolson_scheme
Crank-Nicolson dynamic solver.
subroutine, public solver_dynamic_degree_set(SOLVER, DEGREE, ERR, ERROR,)
Sets/changes the degree of the polynomial used to interpolate time for a dynamic solver.
This module handles all equations matrix and rhs routines.
integer(intg), parameter, public solver_dynamic_first_order
Dynamic solver has first order terms.
subroutine, public solver_type_set(SOLVER, SOLVE_TYPE, ERR, ERROR,)
Sets/changes the type for a solver.
Contains information on an equations set.
This module handles all equations routines.
This module contains all string manipulation and transformation routines.
subroutine, public solvers_create_start(CONTROL_LOOP, SOLVERS, ERR, ERROR,)
Start the creation of a solvers for the control loop.
Contains information on the solvers to be used in a control loop.
subroutine, public fsi_problemspecificationset(problem, problemSpecification, err, error,)
Sets the problem specification for a finite elasticity Navier-Stokes equation type.
subroutine, public fsi_control_loop_post_loop(ControlLoop, Err, Error,)
Runs after each control loop iteration. Updates interface and fluid geometric fields and exports fiel...
Contains information for the interface condition data.
subroutine, public fsi_equations_set_solution_method_set(EQUATIONS_SET, SOLUTION_METHOD, Err, Error,)
Sets/changes the solution method for a finite elasticity navier stokes equation type of a multi physi...
subroutine, public solver_dynamic_linearity_type_set(SOLVER, LINEARITY_TYPE, ERR, ERROR,)
Sets/changes the linearity type for the dynamic solver.
subroutine, public fsi_equationssetspecificationset(equationsSet, specification, Err, Error,)
Sets/changes the equation subtype for a finite elasticity navier stokes equation type of a multi phys...
subroutine, public fsi_post_solve(ControlLoop, Solver, Err, Error,)
Sets up the finite elasticity navier stokes problem post solve.
integer(intg), parameter solver_equations_static
Solver equations are static.
subroutine, public solver_equations_sparsity_type_set(SOLVER_EQUATIONS, SPARSITY_TYPE, ERR, ERROR,)
Sets/changes the sparsity type for solver equations.
subroutine, public solvers_solver_get(SOLVERS, SOLVER_INDEX, SOLVER, ERR, ERROR,)
Returns a pointer to the specified solver in the list of solvers.
Contains information for a field defined on a region.
integer(intg), parameter equations_set_fluid_mechanics_class
integer(intg), parameter solver_equations_linear
Solver equations are linear.
Contains information on a control loop.
subroutine, public fsi_problem_setup(PROBLEM, PROBLEM_SETUP, Err, Error,)
Sets up the finite elasticity navier stokes equations problem.
This module provides an iso_varying_string module, conformant to the API specified in ISO/IEC 1539-2:...
subroutine, public solver_equations_create_finish(SOLVER_EQUATIONS, ERR, ERROR,)
Finishes the process of creating solver equations.
integer(intg), parameter, public solver_sparse_matrices
Use sparse solver matrices.
subroutine, public solver_equations_create_start(SOLVER, SOLVER_EQUATIONS, ERR, ERROR,)
Starts the process of creating solver equations.
integer(intg), parameter, public solver_dynamic_type
A dynamic solver.
integer(intg), parameter problem_setup_solvers_type
Solver setup for a problem.
integer(intg), parameter equations_set_ale_navier_stokes_subtype
This module contains all program wide constants.
integer(intg), parameter solver_equations_nonlinear
Solver equations are nonlinear.
subroutine, public solver_library_type_set(SOLVER, SOLVER_LIBRARY_TYPE, ERR, ERROR,)
Sets/changes the type of library type to use for the solver.
integer(intg), parameter equations_set_finite_elasticity_navier_stokes_ale_subtype
subroutine, public fsi_pre_solve(ControlLoop, Solver, Err, Error,)
Sets up the finite elasticity navier stokes problem pre-solve.
integer(intg), parameter problem_setup_initial_type
Initial setup for a problem.
subroutine, public fsi_equations_set_setup(EQUATIONS_SET, EQUATIONS_SET_SETUP, Err, Error,)
Sets up the finite elasticity navier stokes equation.
subroutine, public solver_equations_linearity_type_set(SOLVER_EQUATIONS, LINEARITY_TYPE, ERR, ERROR,)
Sets/changes the linearity type for solver equations.
subroutine, public exits(NAME)
Records the exit out of the named procedure.
recursive subroutine, public control_loop_solvers_get(CONTROL_LOOP, SOLVERS, ERR, ERROR,)
Returns a pointer to the solvers for a control loop.
integer(intg), parameter equations_set_mooney_rivlin_subtype
This module contains all type definitions in order to avoid cyclic module references.
integer(intg), parameter equations_set_elasticity_class
This module contains all the low-level base routines e.g., all debug, control, and low-level communic...
integer(intg), parameter, public general_output_type
General output type.
Returns the specified control loop as indexed by the control loop identifier from the control loop ro...
subroutine, public control_loop_type_set(CONTROL_LOOP, LOOP_TYPE, ERR, ERROR,)
Sets/changes the control loop type.
integer(intg), parameter problem_multi_physics_class
integer(intg), parameter problem_setup_finish_action
Finish setup action.
This module handles all equations mapping routines.
Contains information about the solver equations for a solver.
integer(intg), parameter equations_set_gfv_solution_method
Grid-based Finite Volume solution method.
Contains information for a problem.
This module handles all Navier-Stokes fluid routines.
This module handles all solver routines.
integer(intg), parameter equations_set_finite_elasticity_type
integer(intg), parameter problem_finite_elasticity_navier_stokes_ale_subtype
integer(intg), parameter problem_finite_elasticity_navier_stokes_type
Implements lists of Field IO operation.
subroutine, public control_loop_create_start(PROBLEM, CONTROL_LOOP, ERR, ERROR,)
Start the process of creating a control loop for a problem.
integer(intg), parameter problem_setup_solver_equations_type
Solver equations setup for a problem.
Contains information on the solver mapping between the global equation sets and the solver matrices...
subroutine, public control_loop_output_type_set(CONTROL_LOOP, OUTPUT_TYPE, ERR, ERROR,)
Sets/changes the output type for a control loop.
subroutine, public solver_dynamic_scheme_set(SOLVER, SCHEME, ERR, ERROR,)
Sets/changes the scheme for a dynamic solver.
Contains information on the nodes defined on a region.
subroutine, public fsi_control_loop_pre_loop(ControlLoop, Err, Error,)
Runs before each control loop iteration.
integer(intg), parameter, public solver_dynamic_nonlinear
Dynamic solver has nonlinear terms.
integer(intg), parameter equations_set_fd_solution_method
Finite Difference solution method.
Contains information on the setup information for an equations set.
integer(intg), parameter problem_setup_start_action
Start setup action.
subroutine, public solver_equations_time_dependence_type_set(SOLVER_EQUATIONS, TIME_DEPENDENCE_TYPE, ERR, ERROR,)
Sets/changes the time dependence type for solver equations.
integer(intg), parameter equations_set_compressible_finite_elasticity_subtype
This module handles all control loop routines.
integer(intg), parameter, public solver_cmiss_library
CMISS (internal) solver library.
subroutine, public errors(NAME, ERR, ERROR)
Records the exiting error of the subroutine.
This module defines all constants shared across equations set routines.
Contains information for the interface data.
integer(intg), parameter equations_set_bem_solution_method
Boundary Element Method solution method.
subroutine, public solver_solver_equations_get(SOLVER, SOLVER_EQUATIONS, ERR, ERROR,)
Returns a pointer to the solver equations for a solver.
integer(intg), parameter equations_set_fv_solution_method
Finite Volume solution method.
integer(intg), parameter, public solver_dynamic_first_degree
Dynamic solver uses a first degree polynomial for time interpolation.
recursive subroutine, public control_loop_create_finish(CONTROL_LOOP, ERR, ERROR,)
Finish the process of creating a control loop.
Flags an error condition.
This module handles all routines pertaining to finite elasticity coupled with navier stokes for fsi p...
integer(intg), parameter, public solver_linear_type
A linear solver.
This module handles all finite elasticity routines.
subroutine, public field_io_elements_export(FIELDS, FILE_NAME, METHOD, ERR, ERROR,)
Export elemental information into multiple files.
subroutine, public fsi_finite_element_calculate(EQUATIONS_SET, ELEMENT_NUMBER, Err, Error,)
Calculates the element stiffness matrices and RHS for a finite elasticity navier stokes equation fini...
This module contains all kind definitions.
subroutine, public field_io_nodes_export(FIELDS, FILE_NAME, METHOD, ERR, ERROR,)
Export nodal information.
integer(intg), parameter equations_set_navier_stokes_equation_type