FEDEM Solver  R8.0
Source code of the dynamics solver
Functions/Subroutines | Variables
solvermodule Module Reference

Module with model containers and top lever driver subroutines. More...

Functions/Subroutines

subroutine, public initialize (chfsi, stateData, ndat, xinp, nxinp, ierr)
 Initializes the simulation model. More...
 
subroutine, public softrestart (stateData, ndat, writeToRDB, ierr)
 Restarts a running simulation model from the provided state. More...
 
subroutine, public solvestep (iop, finalStep, finished, ierr)
 Solves for the next time step or iteration. More...
 
subroutine, public solverampup (iop, finished, ierr)
 Solves the ramp-up stage, if any. More...
 
subroutine, public solvemodes (nModes, eVal, eVec, dofOrder, useLaPack, ierr)
 Solves the eigenvalue system at current configuration. More...
 
subroutine, public finalize (ierr)
 Terminates the simulation and close the result database. More...
 
subroutine, public closeall (ierr, abortsub)
 Closes the result database files and any external modules used. More...
 
subroutine writeclosure (lcon, ierr)
 Outputs file and profiling information on program termination. More...
 
subroutine deallocateall ()
 Deallocates all dynamically allocated data. More...
 
logical function, public solvedynamic ()
 Returns whether the next step should be solved dynamically or not. More...
 
subroutine, public gettime (time, nextStep, ierr)
 Returns the simulation time of the current or next time step. More...
 
subroutine, public setnewtime (nextTime, ierr)
 Sets the time increment size to be used for next time step. More...
 
subroutine, public getengine (value, ierr, userId, tag, x)
 Returns the current value of the specified engine. More...
 
integer function, public getengineid (tag)
 Returns the user Id of the (first) engine with the specified tag. More...
 
subroutine, public getsystemmatrix (Nmat, iopM, ierr)
 Returns a system matrix as a full matrix. More...
 
subroutine, public getelementmatrix (Emat, baseId, iopM, ierr)
 Returns an element matrix for the specified superelement. More...
 
subroutine, public getrhsvector (Rvec, iopV, ierr)
 Returns the current right-hand-side vector of the linearized system. More...
 
subroutine, public setrhsvector (Rvec, addTo, keepDuringIterations, ierr)
 Sets/Updates current right-hand-side vector of linearized system. More...
 
subroutine, public systemsize (ndim, expanded)
 Returns the dimension of the linearized system. More...
 
subroutine, public solverparameters (dt, alpha, beta, gamma)
 Returns parameters for the Newmark/HHT-algorithm. More...
 
subroutine, public objectequations (baseId, meqn, ndof)
 Returns the equation numbers for the specified object. More...
 
subroutine, public objectstatevar (baseId, vars, nVar)
 Returns the current state variables for the specified object. More...
 
subroutine, public solvelineqsystem (rhs, nrhs, ierr)
 Solves current linear equation system for a set of right-hand-sides. More...
 
subroutine, public statevectorsize (transOnly, ndat)
 Returns the dimension of the state vector. More...
 
subroutine, public savestate (transOnly, stateData, ndat, ierr)
 Saves current state to an in-core array. More...
 
subroutine, public partstatevectorsize (iopS, bid, ndat)
 Returns the state vector dimension for a FE part. More...
 
subroutine, public savepartstate (iopS, bid, data, ndat, ierr)
 Saves the deformation/stress state of a FE part to an in-core array. More...
 
subroutine, public straingagessize (ndat)
 Returns the dimension of the strain gages vector. More...
 
subroutine, public saveinitgagestrains (iopS, data, ndat, ierr)
 Save/restores the initial gage strain to/from an in-core array. More...
 
subroutine, public computegagestrains (disp, gageIds, nDisp, nGage, eps, ierr)
 Computes strain tensor at gage positions for input displacements. More...
 
subroutine, public computebeamforces (disp, beamIds, nDofs, nBeams, forces, ierr)
 Computes beam sectional forces for input displacements. More...
 
subroutine, public computerelativedistance (disp, Ids, relDis, nDofs, nIds, ierr)
 Computes relative distance between 2 triads for input displacements. More...
 
subroutine, public computeresponsevars (disp, Ids, resp, nDofs, nIds, ierr)
 Computes response variables for input displacements. More...
 
subroutine, public getjointspringstiffness (sprCoeff, bid, ierr)
 Gets joint spring stiffness coefficients. More...
 

Variables

type(samtype), save, public sam
 Data for managing system matrix assembly. More...
 
type(modestype), save modes
 Data for eigenmodes. More...
 
type(systemtype), save, public sys
 System level model data. More...
 
type(controltype), save ctrl
 Control system data. More...
 
type(mechanismtype), save, public mech
 Mechanism objects of the model. More...
 
type(sysmatrixtype), save amat
 System matrix for external communication. More...
 
real(dp), dimension(:), allocatable, save extrhs
 Externally set load vector. More...
 
character(len=lfnam_p), save chmodel
 Model file name. More...
 
character(len=lfnam_p), dimension(5), save chnames
 Result database file names. More...
 
character(len=lfnam_p), save yamlfile
 Mode shape file name. More...
 
character(len=:), allocatable, save frsnames
 Recovery result files. More...
 

Detailed Description

Module with model containers and top lever driver subroutines.

This module contains all model data and the top-level subroutines of the FEDEM Dynamics solver. It replaces the old top-level driver subroutine of version R7.2 and earlier. The new driver is divided into three main parts:

In addition, this module contains several utility subroutines for external manipulation of the solution process through the solver API.

Function/Subroutine Documentation

◆ closeall()

subroutine, public solvermodule::closeall ( integer, intent(inout)  ierr,
character(len=*), intent(in), optional  abortsub 
)

Closes the result database files and any external modules used.

Parameters
ierrError flag
[in]abortsubIf present, gives the name of the aborting subroutine
Author
Knut Morten Okstad
Date
30 Nov 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ computebeamforces()

subroutine, public solvermodule::computebeamforces ( real(dp), dimension(*), intent(in)  disp,
integer, dimension(*), intent(in)  beamIds,
integer, intent(in)  nDofs,
integer, intent(in)  nBeams,
real(dp), dimension(*), intent(out)  forces,
integer, intent(out)  ierr 
)

Computes beam sectional forces for input displacements.

Parameters
[in]dispDisplacement vector
[in]beamIdsBase IDs of the beams/triads to calculate for
[in]nDofsLength of displacement vector
[in]nBeamsLength of beam identifier array
[out]forcesBeam sectional forces
[out]ierrError flag

The array beamIds contains triplets of (beamId, triadId, dof), where triadId is the base ID of one of the two end triads of the beam element, and dof is a zero-based force component index (-1 implies all). This subroutine always returns all force components, i.e., the dof index is not used here.

Author
Guenter Glanzer
Date
20 Sep 2020
Here is the call graph for this function:
Here is the caller graph for this function:

◆ computegagestrains()

subroutine, public solvermodule::computegagestrains ( real(dp), dimension(*), intent(in)  disp,
integer, dimension(*), intent(in)  gageIds,
integer, intent(in)  nDisp,
integer, intent(in)  nGage,
real(dp), dimension(*), intent(out)  eps,
integer, intent(out)  ierr 
)

Computes strain tensor at gage positions for input displacements.

Parameters
[in]dispDisplacement vector
[in]gageIdsBase IDs of strain gages to calculate for
[in]nDispLength of displacement vector (if 0, use internal state)
[in]nGageLength of gageIds array (number of strain gages)
[out]epsStrain tensors at gage positions
[out]ierrError flag
Author
Guenter Glanzer
Date
9 Sep 2020
Here is the call graph for this function:
Here is the caller graph for this function:

◆ computerelativedistance()

subroutine, public solvermodule::computerelativedistance ( real(dp), dimension(*), intent(in)  disp,
integer, dimension(*), intent(in)  Ids,
real(dp), dimension(*), intent(out)  relDis,
integer, intent(in)  nDofs,
integer, intent(in)  nIds,
integer, intent(out)  ierr 
)

Computes relative distance between 2 triads for input displacements.

Parameters
[in]dispDisplacement vector
[in]IdsUser IDs of the functions with relative sensors
[out]relDisRelative distances
[in]nDofsLength of displacement vector
[in]nIdsLength of function identifier array
[out]ierrError flag
Author
Guenter Glanzer
Date
22 Oct 2020
Here is the call graph for this function:
Here is the caller graph for this function:

◆ computeresponsevars()

subroutine, public solvermodule::computeresponsevars ( real(dp), dimension(*), intent(in)  disp,
integer, dimension(*), intent(in)  Ids,
real(dp), dimension(*), intent(out)  resp,
integer, intent(in)  nDofs,
integer, intent(in)  nIds,
integer, intent(out)  ierr 
)

Computes response variables for input displacements.

Parameters
[in]dispDisplacement vector
[in]IdsUser IDs of the functions with response quantities
[out]respResponse variable values
[in]nDofsLength of displacement vector
[in]nIdsLength of function identifier array
[out]ierrError flag
Author
Guenter Glanzer
Date
6 June 2023
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deallocateall()

subroutine solvermodule::deallocateall
private

Deallocates all dynamically allocated data.

Author
Knut Morten Okstad
Date
23 Jan 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalize()

subroutine, public solvermodule::finalize ( integer, intent(out)  ierr)

Terminates the simulation and close the result database.

Parameters
[out]ierrError flag

This subroutine is invoked after the time integration has completed successfully, and there are no more calculations to be performed. Its main task is to close down the result database, but it also performs the automatic curve export when that has been requested.

Author
Knut Morten Okstad
Date
30 Nov 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getelementmatrix()

subroutine, public solvermodule::getelementmatrix ( real(dp), dimension(*), intent(out)  Emat,
integer, intent(in)  baseId,
integer, intent(in)  iopM,
integer, intent(out)  ierr 
)

Returns an element matrix for the specified superelement.

Parameters
[out]EmatThe element matrix
[in]baseIdBase ID of the superelement to return the matrix for
[in]iopMOption telling which matrix to return (see below)
[out]ierrError flag

The value of iopM is iterpreted as follows:

  • = 1 : Return the element stiffness matrix
  • = 2 : Return the element mass matrix
  • = 3 : Return the element damping matrix
Author
Knut Morten Okstad
Date
16 Feb 2018
Here is the caller graph for this function:

◆ getengine()

subroutine, public solvermodule::getengine ( real(dp), intent(out)  value,
integer, intent(inout)  ierr,
integer, intent(in)  userId,
character(len=*), intent(in), optional  tag,
real(dp), intent(in), optional  x 
)

Returns the current value of the specified engine.

Parameters
[out]valueThe function value
ierrError flag
[in]userIdUser ID of the engine to evaluate
[in]xOptional function argument value
[in]tagOptional tag to identify the engine with

If an engine with userId is found, the value of that engine is returned. Otherwise, the (first) engine with a matching tag is used. If no engine with the given identification is found, zero is returned. The error flag is decrementend on failure, othwerwise it is not touched.

Author
Knut Morten Okstad
Date
5 Dec 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getengineid()

integer function, public solvermodule::getengineid ( character(len=*), intent(in)  tag)

Returns the user Id of the (first) engine with the specified tag.

Parameters
[in]tagTag to identify the engine with
Returns
User ID of tagged function, or -1 if not found
Author
Knut Morten Okstad
Date
25 Jun 2022
Here is the caller graph for this function:

◆ getjointspringstiffness()

subroutine, public solvermodule::getjointspringstiffness ( real(dp), dimension(*), intent(out)  sprCoeff,
integer, intent(in)  bid,
integer, intent(out)  ierr 
)

Gets joint spring stiffness coefficients.

Parameters
[out]sprCoeffSpring stiffness coefficient
[in]bidBase ID for spring
[out]ierrError flag
Author
Guenter Glanzer
Date
20 June 2023
Here is the caller graph for this function:

◆ getrhsvector()

subroutine, public solvermodule::getrhsvector ( real(dp), dimension(sam%neq), intent(out)  Rvec,
integer, intent(in)  iopV,
integer, intent(out)  ierr 
)

Returns the current right-hand-side vector of the linearized system.

Parameters
[out]RvecThe right-hand-side vector
[in]iopVIf equal to 1, return current external load vector instead
[out]ierrError flag
Author
Knut Morten Okstad
Date
30 Mar 2017
Here is the caller graph for this function:

◆ getsystemmatrix()

subroutine, public solvermodule::getsystemmatrix ( real(dp), dimension(sam%neq,sam%neq), intent(out)  Nmat,
integer, intent(in)  iopM,
integer, intent(out)  ierr 
)

Returns a system matrix as a full matrix.

Parameters
[out]NmatThe system matrix
[in]iopMOption telling which matrix to return (see below)
[out]ierrError flag

The value of iopM is iterpreted as follows:

  • = 1 : Return the system stiffness matrix
  • = 2 : Return the system mass matrix
  • = 3 : Return the system damping matrix
  • > 3 : Return the current system Newton matrix
Author
Knut Morten Okstad
Date
30 Mar 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gettime()

subroutine, public solvermodule::gettime ( real(dp), intent(out)  time,
logical, intent(in)  nextStep,
integer, intent(inout)  ierr 
)

Returns the simulation time of the current or next time step.

Parameters
[out]timeThe simulation time returned
[in]nextStepIf .true., return for next step, otherwise current step
ierrError flag

This subroutine does not always work with nextStep = .true. if cut-back is performed, as the time step size then might be decreased.

Author
Knut Morten Okstad
Date
5 Dec 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialize()

subroutine, public solvermodule::initialize ( character(len=*), intent(in), optional  chfsi,
real(dp), dimension(*), intent(in), optional  stateData,
integer, intent(in), optional  ndat,
real(dp), dimension(*), intent(in), optional  xinp,
integer, intent(inout), optional  nxinp,
integer, intent(out)  ierr 
)

Initializes the simulation model.

Parameters
[in]chfsiText string containing the model definition
[in]stateDataCurrent state to restart from
[in]ndatLength of the stateData array
[in]xinpExternal function values for initial state
nxinpNumber of external function values to extract/extracted
[out]ierrError flag

This subroutine reads and sets up a new simulation model, defines the initial configuration and initializes the result database. It also handles restart, both from file and from a solution state provided as argument.

Author
Knut Morten Okstad
Date
30 Nov 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ objectequations()

subroutine, public solvermodule::objectequations ( integer, intent(in)  baseId,
integer, dimension(*), intent(out)  meqn,
integer, intent(out)  ndof 
)

Returns the equation numbers for the specified object.

Parameters
[in]baseIdBase ID of the object to get equation numbers for
[out]meqnArray of equation numbers
[out]ndofNumber of degrees of freedom for the object in question
Author
Knut Morten Okstad
Date
30 Mar 2017
Here is the caller graph for this function:

◆ objectstatevar()

subroutine, public solvermodule::objectstatevar ( integer, intent(in)  baseId,
real(dp), dimension(*), intent(out)  vars,
integer, intent(out)  nVar 
)

Returns the current state variables for the specified object.

Parameters
[in]baseIdBase ID of the object to get equation numbers for
[out]varsArray of state variables (positions only for Triads)
[out]nVarNumber of state variables
Author
Knut Morten Okstad
Date
20 Mar 2020
Here is the caller graph for this function:

◆ partstatevectorsize()

subroutine, public solvermodule::partstatevectorsize ( integer, intent(in)  iopS,
integer, intent(in)  bid,
integer, intent(out)  ndat 
)

Returns the state vector dimension for a FE part.

Parameters
[in]iopSOption telling for which state, 1: deformation, 2: stresses
[in]bidBase ID if the FE part to consider
[in]ndatRequired length of the state array

The deformation state size equals 4 + 3*(number of nodes) whereas the stress state size euals 4 + number of nodes. In both cases the first 4 entries are the step number, time, step size and baseID, respectiviely.

Author
Runar Heggelien Refsnaes
Date
29 Oct 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ saveinitgagestrains()

subroutine, public solvermodule::saveinitgagestrains ( integer, intent(in)  iopS,
real(dp), dimension(*), intent(inout)  data,
integer, intent(in)  ndat,
integer, intent(out)  ierr 
)

Save/restores the initial gage strain to/from an in-core array.

Parameters
[in]iopSIf equal to 1 save to data array, otherwise restore from it
dataArray with strain gage values
[in]ndatLength of the data array
[out]ierrError flag
Author
Knut Morten Okstad
Date
18 Dec 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ savepartstate()

subroutine, public solvermodule::savepartstate ( integer, intent(in)  iopS,
integer, intent(in)  bid,
real(dp), dimension(*), intent(out)  data,
integer, intent(in)  ndat,
integer, intent(out)  ierr 
)

Saves the deformation/stress state of a FE part to an in-core array.

Parameters
[in]iopSOption telling what to save, 1: deformation, 2: stress
[in]bidBase ID if the FE part to consider
[out]dataDeformation/stress state array
[in]ndatLength of the data array
[out]ierrError flag

If iopS equals 1, the von Mises stress value at each node of the part is saved, otherwise the 3 deformational displacement components at each node are saved. As the first four values, the current step number, the current time, the time step size, and the base ID of the FE part, respectively, are saved.

Author
Runar Heggelien Refsnaes
Date
31 Oct 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ savestate()

subroutine, public solvermodule::savestate ( logical, intent(in)  transOnly,
real(dp), dimension(*), intent(out)  stateData,
integer, intent(in)  ndat,
integer, intent(out)  ierr 
)

Saves current state to an in-core array.

Parameters
[in]transOnlyIf .true., consider position matrices only
[out]stateDataArray containing all solution dependent variables
[in]ndatLength of the stateData array
[out]ierrError flag
Author
Knut Morten Okstad
Date
16 Jan 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setnewtime()

subroutine, public solvermodule::setnewtime ( real(dp), intent(in)  nextTime,
integer, intent(out)  ierr 
)

Sets the time increment size to be used for next time step.

Parameters
[in]nextTimePhysical time to calculate time increment size from
[out]ierrError flag
Author
Knut Morten Okstad
Date
22 Jun 2022
Here is the caller graph for this function:

◆ setrhsvector()

subroutine, public solvermodule::setrhsvector ( real(dp), dimension(sam%neq), intent(in)  Rvec,
logical, intent(in)  addTo,
logical, intent(in)  keepDuringIterations,
integer, intent(out)  ierr 
)

Sets/Updates current right-hand-side vector of linearized system.

Parameters
RvecThe right-hand-side vector
[in]addToIf .true., the provided vector is added to the current right-hand-side vector in the model, otherwise it replaces the current one
[in]keepDuringIterationsIf .true., the provided vector should be applied in all iterations of current time step
[out]ierrError flag
Author
Knut Morten Okstad
Date
30 Mar 2017
Here is the caller graph for this function:

◆ softrestart()

subroutine, public solvermodule::softrestart ( real(dp), dimension(*), intent(in)  stateData,
integer, intent(in)  ndat,
integer, intent(in)  writeToRDB,
integer, intent(out)  ierr 
)

Restarts a running simulation model from the provided state.

Parameters
[in]stateDataCurrent state to restart from
[in]ndatLength of the stateData array
[in]writeToRDBControl variable for frs-file output
[out]ierrError flag

The value of the argument writeToRDB is interpreted as follows:

  • < 1 : Suppress all frs-output
  • = 1 : Continue with output to same files as the previous run
  • = 2 : Create a new set of files for this restart run
Author
Knut Morten Okstad
Date
6 Feb 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solvedynamic()

logical function, public solvermodule::solvedynamic

Returns whether the next step should be solved dynamically or not.

Author
Guenter Glanzer
Date
10 Jun 2018
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solvelineqsystem()

subroutine, public solvermodule::solvelineqsystem ( real(dp), dimension(:), intent(inout)  rhs,
integer, intent(in)  nrhs,
integer, intent(out)  ierr 
)

Solves current linear equation system for a set of right-hand-sides.

Parameters
rhsThe right-hand-side vectors of the linear equation system
[in]nrhsNumber of right-hand-side vectors
[out]ierrError flag

Assuming the Newton matrix already have been factorized. To be used by the inverse problem algorithm.

Author
Knut Morten Okstad
Date
16 Jan 2017
Here is the caller graph for this function:

◆ solvemodes()

subroutine, public solvermodule::solvemodes ( integer, intent(in)  nModes,
real(dp), dimension(*), intent(out)  eVal,
real(dp), dimension(*), intent(out)  eVec,
logical, intent(in)  dofOrder,
integer, intent(in)  useLaPack,
integer, intent(out)  ierr 
)

Solves the eigenvalue system at current configuration.

Parameters
[in]nModesNumber of eigenmodes to solve for
[out]eValComputed eigenvalues (frequencies in [Hz])
[out]eVecMass-normalized eigenvectors
[in]dofOrderIf .true., expand eigenvectors to DOF-order
[in]useLaPackFlag for using LAPack eigenvalue solvers
[out]ierrError flag
Author
Knut Morten Okstad
Date
26 Oct 2020
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solverampup()

subroutine, public solvermodule::solverampup ( integer, intent(inout)  iop,
logical, intent(out)  finished,
integer, intent(out)  ierr 
)

Solves the ramp-up stage, if any.

Parameters
iopOperation flag telling what to do, should be zero
[out]finishedIf .true., the time integration failed
[out]ierrError flag
Author
Knut Morten Okstad
Date
19 Aug 2022
Here is the call graph for this function:
Here is the caller graph for this function:

◆ solverparameters()

subroutine, public solvermodule::solverparameters ( real(dp), intent(out)  dt,
real(dp), intent(out)  alpha,
real(dp), intent(out)  beta,
real(dp), intent(out)  gamma 
)

Returns parameters for the Newmark/HHT-algorithm.

Parameters
[out]dtTime step size
[out]alphaNumerical damping parameter of the HHT-scheme
[out]betaNewmark time integration parameter
[out]gammaNewmark time integration parameter
Author
Guenter Glanzer
Date
10 Jun 2018
Here is the caller graph for this function:

◆ solvestep()

subroutine, public solvermodule::solvestep ( integer, intent(inout)  iop,
logical, intent(in)  finalStep,
logical, intent(out)  finished,
integer, intent(out)  ierr 
)

Solves for the next time step or iteration.

Parameters
iopControl variable defining what to do (see below)
[in]finalStepIf .true., this is the final step of current time window
[out]finishedIf .true., the end of the simulation has been reached, or the time integration failed
[out]ierrError flag

This subroutine performs different sub-tasks of the solution process, depending on the value of the control variable iop, as follows:

  • = 0 : Normal operation, solve the entire step with iterations
  • = -1 : Establish the linear system of the first (prediction) iteration, and then return to the calling module
  • = -2 : Same as -1, but also factorize the Newton matrix
  • = -3 : Same as 0, but this is the first time step after restart
  • = 1 : We are doing quasi-static iterations (internal mode)
  • = 2 : We are doing Newmark iterations (internal mode)
  • = 3 : Do one quasi-static iteration on the existing linear system
  • = 4 : Do one Newmark iteration on the existing linear system
  • = 5 : As 3, and then continue until equilibrium
  • = 6 : As 4, and then continue until equilibrium
Author
Knut Morten Okstad
Date
30 Nov 2016
Here is the call graph for this function:
Here is the caller graph for this function:

◆ statevectorsize()

subroutine, public solvermodule::statevectorsize ( logical, intent(in)  transOnly,
integer, intent(out)  ndat 
)

Returns the dimension of the state vector.

Parameters
[in]transOnlyIf .true., consider position matrices only
[out]ndatRequired length of the state array
Author
Knut Morten Okstad
Date
16 Jan 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ straingagessize()

subroutine, public solvermodule::straingagessize ( integer, intent(out)  ndat)

Returns the dimension of the strain gages vector.

Parameters
[out]ndatRequired length of the strain gage state array
Author
Knut Morten Okstad
Date
18 Dec 2017
Here is the call graph for this function:
Here is the caller graph for this function:

◆ systemsize()

subroutine, public solvermodule::systemsize ( integer, intent(out)  ndim,
logical, intent(in), optional  expanded 
)

Returns the dimension of the linearized system.

Parameters
[out]ndimNumber of unknowns in the linear equation system
[in]expandedIf .true., return total number of DOFs in the system
Author
Knut Morten Okstad
Date
30 Mar 2017
Here is the caller graph for this function:

◆ writeclosure()

subroutine solvermodule::writeclosure ( integer, intent(in)  lcon,
integer, intent(inout)  ierr 
)

Outputs file and profiling information on program termination.

Parameters
[in]lconFile unit number for console messages
ierrError flag
Author
Knut Morten Okstad
Date
30 Nov 2016
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ amat

type(sysmatrixtype), save solvermodule::amat
private

System matrix for external communication.

◆ chmodel

character(len=lfnam_p), save solvermodule::chmodel
private

Model file name.

◆ chnames

character(len=lfnam_p), dimension(5), save solvermodule::chnames
private

Result database file names.

◆ ctrl

type(controltype), save solvermodule::ctrl
private

Control system data.

◆ extrhs

real(dp), dimension(:), allocatable, save solvermodule::extrhs
private

Externally set load vector.

Assumed constant during Newton iterations.

◆ frsnames

character(len=:), allocatable, save solvermodule::frsnames
private

Recovery result files.

◆ mech

type(mechanismtype), save, public solvermodule::mech

Mechanism objects of the model.

◆ modes

type(modestype), save solvermodule::modes
private

Data for eigenmodes.

◆ sam

type(samtype), save, public solvermodule::sam

Data for managing system matrix assembly.

◆ sys

type(systemtype), save, public solvermodule::sys

System level model data.

◆ yamlfile

character(len=lfnam_p), save solvermodule::yamlfile
private

Mode shape file name.