v0.14.0
Public Member Functions | Public Attributes | List of all members
ConvectiveMassElement::ShellResidualElement Struct Reference

#include <users_modules/basic_finite_elements/src/ConvectiveMassElement.hpp>

Inheritance diagram for ConvectiveMassElement::ShellResidualElement:
[legend]
Collaboration diagram for ConvectiveMassElement::ShellResidualElement:
[legend]

Public Member Functions

 ShellResidualElement (MoFEM::Interface &m_field)
 
MoFEMErrorCode preProcess ()
 Calculate inconsistency between approximation of velocities and velocities calculated from displacements. More...
 

Public Attributes

MoFEM::InterfacemField
 
MatShellCtxshellMatCtx
 pointer to shell matrix More...
 

Detailed Description

Examples
nonlinear_dynamics.cpp.

Definition at line 715 of file ConvectiveMassElement.hpp.

Constructor & Destructor Documentation

◆ ShellResidualElement()

ConvectiveMassElement::ShellResidualElement::ShellResidualElement ( MoFEM::Interface m_field)

Definition at line 2405 of file ConvectiveMassElement.cpp.

2407  : mField(m_field) {}

Member Function Documentation

◆ preProcess()

MoFEMErrorCode ConvectiveMassElement::ShellResidualElement::preProcess ( )

Calculate inconsistency between approximation of velocities and velocities calculated from displacements.

Returns
MoFEMErrorCode

Definition at line 2409 of file ConvectiveMassElement.cpp.

2409  {
2411 
2412  if (ts_ctx != CTX_TSSETIFUNCTION) {
2413  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
2414  "It is used to residual of velocities");
2415  }
2416  if (!shellMatCtx->iNitialized) {
2417  CHKERR shellMatCtx->iNit();
2418  }
2419  // Note velocities calculate from displacements are stroed in shellMatCtx->u
2420  CHKERR VecScatterBegin(shellMatCtx->scatterU, ts_u_t, shellMatCtx->u,
2421  INSERT_VALUES, SCATTER_FORWARD);
2422  CHKERR VecScatterEnd(shellMatCtx->scatterU, ts_u_t, shellMatCtx->u,
2423  INSERT_VALUES, SCATTER_FORWARD);
2424  CHKERR VecScatterBegin(shellMatCtx->scatterV, ts_u, shellMatCtx->v,
2425  INSERT_VALUES, SCATTER_FORWARD);
2426  CHKERR VecScatterEnd(shellMatCtx->scatterV, ts_u, shellMatCtx->v,
2427  INSERT_VALUES, SCATTER_FORWARD);
2428  CHKERR VecAXPY(shellMatCtx->v, -1, shellMatCtx->u);
2429  CHKERR VecScatterBegin(shellMatCtx->scatterV, shellMatCtx->v, ts_F,
2430  ADD_VALUES, SCATTER_REVERSE);
2431  CHKERR VecScatterEnd(shellMatCtx->scatterV, shellMatCtx->v, ts_F, ADD_VALUES,
2432  SCATTER_REVERSE);
2433  // VecView(shellMatCtx->v,PETSC_VIEWER_STDOUT_WORLD);
2434 
2436 }

Member Data Documentation

◆ mField

MoFEM::Interface& ConvectiveMassElement::ShellResidualElement::mField

Definition at line 716 of file ConvectiveMassElement.hpp.

◆ shellMatCtx

MatShellCtx* ConvectiveMassElement::ShellResidualElement::shellMatCtx

pointer to shell matrix

Definition at line 720 of file ConvectiveMassElement.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
ConvectiveMassElement::MatShellCtx::iNit
MoFEMErrorCode iNit()
Definition: ConvectiveMassElement.cpp:2352
ConvectiveMassElement::MatShellCtx::scatterU
VecScatter scatterU
Definition: ConvectiveMassElement.hpp:505
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
ConvectiveMassElement::ShellResidualElement::mField
MoFEM::Interface & mField
Definition: ConvectiveMassElement.hpp:716
ConvectiveMassElement::MatShellCtx::u
Vec u
Definition: ConvectiveMassElement.hpp:513
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ConvectiveMassElement::ShellResidualElement::shellMatCtx
MatShellCtx * shellMatCtx
pointer to shell matrix
Definition: ConvectiveMassElement.hpp:720
ConvectiveMassElement::MatShellCtx::v
Vec v
Definition: ConvectiveMassElement.hpp:513
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
ConvectiveMassElement::MatShellCtx::scatterV
VecScatter scatterV
Definition: ConvectiveMassElement.hpp:505
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
ConvectiveMassElement::MatShellCtx::iNitialized
bool iNitialized
Definition: ConvectiveMassElement.hpp:508