Mandel's consolidation

A sample's dimensions are and , and it is in plane strain (no displacement). It is squashed with constant normal force by impermeable, frictionless plattens on its top and bottom surfaces (at ). Fluid is allowed to leak out from its sides (at ), but all other surfaces are impermeable. This is called Mandel's problem and it is shown graphically in []

The interesting feature of this problem (apart from that it can be solved analytically) is that the porepressure in the sample's center increases for a short time after application of the force. This is because the leakage of the fluid from the sample's sides causes an apparent softening of the material near those sides. This means stress concentrates towards the sample's center which causes an increase in porepressure. Of course, eventually the fluid totally drains from the sample, and the porepressure is zero. As the fluid drains from the sample's sides the plattens move slowly towards each other.

As is common in the literature, this is simulated by considering the quarter-sample, and , with impermeable, roller BCs at and and . Porepressure is fixed at zero on . Porepressure and displacement are initialised to zero. Then the top () is moved downwards with prescribed velocity, so that the total force that is induced by this downwards velocity is fixed. The velocity is worked out by solving Mandel's problem analytically, and the total force is monitored in the simulation to check that it indeed remains constant.

The simulations in the PorousFlow test suite use 10 elements in the direction and 1 in the direction. Four types of simulation are run:

  1. FullSat. This uses the FullySaturated versions of the fluid mass time derivative and the fluid flux. In this case the Biot modulus is kept fixed, so it is expected to agree with the analytical solutions.

[Mesh<<<{"href": "../syntax/Mesh/index.html"}>>>]
  type = GeneratedMesh
  dim = 3
  nx = 10
  ny = 1
  nz = 1
  xmin = 0
  xmax = 1
  ymin = 0
  ymax = 0.1
  zmin = 0
  zmax = 1
[]

[Variables<<<{"href": "../syntax/Variables/index.html"}>>>]
  [./disp_x]
  [../]
  [./disp_y]
  [../]
  [./disp_z]
  [../]
[]

[BCs<<<{"href": "../syntax/BCs/index.html"}>>>]
  [./roller_xmin]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'left'
  [../]
  [./roller_ymin]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'bottom'
  [../]
  [./plane_strain]
    type = DirichletBC<<<{"description": "Imposes the essential boundary condition $u=g$, where $g$ is a constant, controllable value.", "href": "../source/bcs/DirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
    value<<<{"description": "Value of the BC"}>>> = 0
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = 'back front'
  [../]
  [./top_velocity]
    type = FunctionDirichletBC<<<{"description": "Imposes the essential boundary condition $u=g(t,\\vec{x})$, where $g$ is a (possibly) time and space-dependent MOOSE Function.", "href": "../source/bcs/FunctionDirichletBC.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    function<<<{"description": "The forcing function."}>>> = top_velocity
    boundary<<<{"description": "The list of boundary IDs from the mesh where this object applies"}>>> = top
  [../]
[]

[Functions<<<{"href": "../syntax/Functions/index.html"}>>>]
  [./top_velocity]
    type = PiecewiseLinear<<<{"description": "Linearly interpolates between pairs of x-y data", "href": "../source/functions/PiecewiseLinear.html"}>>>
    x<<<{"description": "The abscissa values"}>>> = '0 0.002 0.006   0.014   0.03    0.046   0.062   0.078   0.094   0.11    0.126   0.142   0.158   0.174   0.19 0.206 0.222 0.238 0.254 0.27 0.286 0.302 0.318 0.334 0.35 0.366 0.382 0.398 0.414 0.43 0.446 0.462 0.478 0.494 0.51 0.526 0.542 0.558 0.574 0.59 0.606 0.622 0.638 0.654 0.67 0.686 0.702'
    y<<<{"description": "The ordinate values"}>>> = '-0.041824842    -0.042730269    -0.043412712    -0.04428867     -0.045509181    -0.04645965     -0.047268246 -0.047974749      -0.048597109     -0.0491467  -0.049632388     -0.050061697      -0.050441198     -0.050776675     -0.051073238      -0.0513354 -0.051567152      -0.051772022     -0.051953128 -0.052113227 -0.052254754 -0.052379865 -0.052490464 -0.052588233 -0.052674662 -0.052751065 -0.052818606 -0.052878312 -0.052931093 -0.052977751 -0.053018997 -0.053055459 -0.053087691 -0.053116185 -0.053141373 -0.05316364 -0.053183324 -0.053200724 -0.053216106 -0.053229704 -0.053241725 -0.053252351 -0.053261745 -0.053270049 -0.053277389 -0.053283879 -0.053289615'
  [../]
[]

[Kernels<<<{"href": "../syntax/Kernels/index.html"}>>>]
  [./grad_stress_x]
    type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../source/kernels/StressDivergenceTensors.html"}>>>
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
    displacements<<<{"description": "The string of displacements suitable for the problem statement"}>>> = 'disp_x disp_y disp_z'
    component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 0
  [../]
  [./grad_stress_y]
    type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../source/kernels/StressDivergenceTensors.html"}>>>
    displacements<<<{"description": "The string of displacements suitable for the problem statement"}>>> = 'disp_x disp_y disp_z'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 1
  [../]
  [./grad_stress_z]
    type = StressDivergenceTensors<<<{"description": "Stress divergence kernel for the Cartesian coordinate system", "href": "../source/kernels/StressDivergenceTensors.html"}>>>
    displacements<<<{"description": "The string of displacements suitable for the problem statement"}>>> = 'disp_x disp_y disp_z'
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
    component<<<{"description": "An integer corresponding to the direction the variable this kernel acts in. (0 for x, 1 for y, 2 for z)"}>>> = 2
  [../]
  [poro_x]
    type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
    PorousFlowDictator<<<{"description": "The UserObject that holds the list of PorousFlow variable names."}>>> = dictator
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.6
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_x
    component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 0
  []
  [poro_y]
    type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
    PorousFlowDictator<<<{"description": "The UserObject that holds the list of PorousFlow variable names."}>>> = dictator
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.6
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_y
    component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 1
  []
  [poro_z]
    type = PorousFlowEffectiveStressCoupling<<<{"description": "Implements the weak form of the expression biot_coefficient * grad(effective fluid pressure)", "href": "../source/kernels/PorousFlowEffectiveStressCoupling.html"}>>>
    PorousFlowDictator<<<{"description": "The UserObject that holds the list of PorousFlow variable names."}>>> = dictator
    biot_coefficient<<<{"description": "Biot coefficient"}>>> = 0.6
    component<<<{"description": "The component (0 for x, 1 for y and 2 for z) of grad(P)"}>>> = 2
    variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = disp_z
  []
[]

[AuxVariables<<<{"href": "../syntax/AuxVariables/index.html"}>>>]
  [./porepressure]
  initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>>=0
  [../]
  [./stress_yy]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  [../]
  [./tot_force]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  [../]
[]

[AuxKernels<<<{"href": "../syntax/AuxKernels/index.html"}>>>]
  [./stress_yy]
    type = RankTwoAux<<<{"description": "Access a component of a RankTwoTensor", "href": "../source/auxkernels/RankTwoAux.html"}>>>
    rank_two_tensor<<<{"description": "The rank two material tensor name"}>>> = stress
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = stress_yy
    index_i<<<{"description": "The index i of ij for the tensor to output (0, 1, 2)"}>>> = 1
    index_j<<<{"description": "The index j of ij for the tensor to output (0, 1, 2)"}>>> = 1
  [../]
  [./tot_force]
    type = ParsedAux<<<{"description": "Sets a field variable value to the evaluation of a parsed expression.", "href": "../source/auxkernels/ParsedAux.html"}>>>
    args = 'stress_yy porepressure'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = timestep_end
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = tot_force
    function = '-stress_yy+0.6*porepressure'
  [../]
[]

[UserObjects<<<{"href": "../syntax/UserObjects/index.html"}>>>]
  [dictator]
    type = PorousFlowDictator<<<{"description": "Holds information on the PorousFlow variable names", "href": "../source/userobjects/PorousFlowDictator.html"}>>>
    porous_flow_vars<<<{"description": "List of primary variables that are used in the PorousFlow simulation.  Jacobian entries involving derivatives wrt these variables will be computed.  In single-phase models you will just have one (eg 'pressure'), in two-phase models you will have two (eg 'p_water p_gas', or 'p_water s_water'), etc."}>>> = 'disp_x disp_y disp_z'
    number_fluid_phases<<<{"description": "The number of fluid phases in the simulation"}>>> = 1
    number_fluid_components<<<{"description": "The number of fluid components in the simulation"}>>> = 1
  []
[]

[Materials<<<{"href": "../syntax/Materials/index.html"}>>>]
  [./elasticity_tensor]
    type = ComputeElasticityTensor<<<{"description": "Compute an elasticity tensor.", "href": "../source/materials/ComputeElasticityTensor.html"}>>>
    C_ijkl<<<{"description": "Stiffness tensor for material"}>>> = '0.5 0.75'
    # bulk modulus is lambda + 2*mu/3 = 0.5 + 2*0.75/3 = 1
    fill_method<<<{"description": "The fill method"}>>> = symmetric_isotropic
  [../]
  [./strain]
    type = ComputeSmallStrain<<<{"description": "Compute a small strain.", "href": "../source/materials/ComputeSmallStrain.html"}>>>
    displacements<<<{"description": "The displacements appropriate for the simulation geometry and coordinate system"}>>> = 'disp_x disp_y disp_z'
  [../]
  [./stress]
    type = ComputeLinearElasticStress<<<{"description": "Compute stress using elasticity for small strains", "href": "../source/materials/ComputeLinearElasticStress.html"}>>>
  [../]

  [vol_strain]
    type = PorousFlowEffectiveFluidPressure<<<{"description": "This Material calculates an effective fluid pressure: effective_stress = total_stress + biot_coeff*effective_fluid_pressure.  The effective_fluid_pressure = sum_{phases}(S_phase * P_phase)", "href": "../source/materials/PorousFlowEffectiveFluidPressure.html"}>>>
    PorousFlowDictator<<<{"description": "The UserObject that holds the list of PorousFlow variable names"}>>> = dictator
  []
  [ppss]
    PorousFlowDictator<<<{"description": "The UserObject that holds the list of PorousFlow variable names"}>>> = dictator
    type = PorousFlow1PhaseFullySaturated<<<{"description": "This Material is used for the fully saturated single-phase situation where porepressure is the primary variable", "href": "../source/materials/PorousFlow1PhaseFullySaturated.html"}>>>
    porepressure<<<{"description": "Variable that represents the porepressure of the single phase"}>>> = porepressure
  []
[]

[Executioner<<<{"href": "../syntax/Executioner/index.html"}>>>]
  type = Transient
  solve_type = Newton
  start_time = 0
  end_time = 0.7

  petsc_options = '-ksp_snes_ew'
  petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
  petsc_options_value = 'lu       superlu_dist'

  #fixed_point_algorithm = 'secant'

  accept_on_max_fixed_point_iteration = true
  fixed_point_max_its = 2
  fixed_point_rel_tol = 1e-20

  relaxation_factor = 0.5
  relaxed_variables = 'disp_x disp_y disp_z'

  [TimeStepper<<<{"href": "../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = PostprocessorDT
    postprocessor = dt
    dt = 0.001
  []
[]

[Postprocessors<<<{"href": "../syntax/Postprocessors/index.html"}>>>]
  [./p0]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.0 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p1]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.1 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p2]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.2 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p3]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.3 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p4]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.4 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p5]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.5 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p6]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.6 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p7]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.7 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p8]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.8 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p9]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '0.9 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./p99]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = porepressure
  [../]
  [./xdisp]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0.1 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_x
  [../]
  [./ydisp]
    type = PointValue<<<{"description": "Compute the value of a variable at a specified location", "href": "../source/postprocessors/PointValue.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
    point<<<{"description": "The physical point where the solution will be evaluated."}>>> = '1 0.1 0'
    variable<<<{"description": "The name of the variable that this postprocessor operates on."}>>> = disp_y
  [../]
  [./total_downwards_force]
     type = ElementAverageValue<<<{"description": "Computes the volumetric average of a variable", "href": "../source/postprocessors/ElementAverageValue.html"}>>>
     outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = csv
     variable<<<{"description": "The name of the variable that this object operates on"}>>> = tot_force
  [../]

  [./dt]
    type = FunctionValuePostprocessor<<<{"description": "Computes the value of a supplied function at a single point (scalable)", "href": "../source/postprocessors/FunctionValuePostprocessor.html"}>>>
    outputs<<<{"description": "Vector of output names where you would like to restrict the output of variables(s) associated with this object"}>>> = console
    function<<<{"description": "The function which supplies the postprocessor value."}>>> = if(0.15*t<0.01,0.15*t,0.01)
  [../]
[]

[Outputs<<<{"href": "../syntax/Outputs/index.html"}>>>]
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = false
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = false
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  file_base<<<{"description": "Common file base name to be utilized with all output objects"}>>> = mandel_master
  [csv]
    interval = 3
    type = CSV<<<{"description": "Output for postprocessors, vector postprocessors, and scalar variables using comma seperated values (CSV).", "href": "../source/outputs/CSV.html"}>>>
  []
[]

[MultiApps<<<{"href": "../syntax/MultiApps/index.html"}>>>]
  [./sub]
    type = TransientMultiApp<<<{"description": "MultiApp for performing coupled simulations with the parent and sub-application both progressing in time.", "href": "../source/multiapps/TransientMultiApp.html"}>>>
    app_type<<<{"description": "The type of application to build (applications not registered can be loaded with dynamic libraries. Parent application type will be used if not provided."}>>> =FalconApp
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = TIMESTEP_END
    input_files<<<{"description": "The input file for each App.  If this parameter only contains one input file it will be used for all of the Apps.  When using 'positions_from_file' it is also admissable to provide one input_file per file."}>>> = mandel_sub.i
  [../]
[]

[Transfers<<<{"href": "../syntax/Transfers/index.html"}>>>]
  [./disp_x_to_sub]
    type = MultiAppCopyTransfer<<<{"description": "Copies variables (nonlinear and auxiliary) between multiapps that have identical meshes.", "href": "../source/transfers/MultiAppCopyTransfer.html"}>>>
    direction<<<{"description": "Whether this Transfer will be 'to' or 'from' a MultiApp, or bidirectional, by providing both FROM_MULTIAPP and TO_MULTIAPP."}>>> = to_multiapp
    multi_app<<<{"description": "The name of the MultiApp to transfer data with"}>>> = sub
    source_variable<<<{"description": "The variable to transfer from."}>>> = disp_x
    variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = disp_x
  [../]
  [./disp_y_to_sub]
    type = MultiAppCopyTransfer<<<{"description": "Copies variables (nonlinear and auxiliary) between multiapps that have identical meshes.", "href": "../source/transfers/MultiAppCopyTransfer.html"}>>>
    direction<<<{"description": "Whether this Transfer will be 'to' or 'from' a MultiApp, or bidirectional, by providing both FROM_MULTIAPP and TO_MULTIAPP."}>>> = to_multiapp
    multi_app<<<{"description": "The name of the MultiApp to transfer data with"}>>> = sub
    source_variable<<<{"description": "The variable to transfer from."}>>> = disp_y
    variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = disp_y
  [../]
  [./disp_z_to_sub]
    type = MultiAppCopyTransfer<<<{"description": "Copies variables (nonlinear and auxiliary) between multiapps that have identical meshes.", "href": "../source/transfers/MultiAppCopyTransfer.html"}>>>
    direction<<<{"description": "Whether this Transfer will be 'to' or 'from' a MultiApp, or bidirectional, by providing both FROM_MULTIAPP and TO_MULTIAPP."}>>> = to_multiapp
    multi_app<<<{"description": "The name of the MultiApp to transfer data with"}>>> = sub
    source_variable<<<{"description": "The variable to transfer from."}>>> = disp_z
    variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = disp_z
  [../]
  [./porepressure_from_sub]
    type = MultiAppCopyTransfer<<<{"description": "Copies variables (nonlinear and auxiliary) between multiapps that have identical meshes.", "href": "../source/transfers/MultiAppCopyTransfer.html"}>>>
    direction<<<{"description": "Whether this Transfer will be 'to' or 'from' a MultiApp, or bidirectional, by providing both FROM_MULTIAPP and TO_MULTIAPP."}>>> = from_multiapp
    multi_app<<<{"description": "The name of the MultiApp to transfer data with"}>>> = sub
    source_variable<<<{"description": "The variable to transfer from."}>>> = porepressure
    variable<<<{"description": "The auxiliary variable to store the transferred values in."}>>> = porepressure
  [../]
[]
(examples/fixedStressMultiApp/mandel_master.i)

Figure 1: The horiztonal and vertical displacement and the pore pressure at different locations of the platten as a function of time.