v0.14.0
Public Member Functions | Public Attributes | List of all members
SimpleArcLengthControl Struct Reference

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

Inheritance diagram for SimpleArcLengthControl:
[legend]
Collaboration diagram for SimpleArcLengthControl:
[legend]

Public Member Functions

 SimpleArcLengthControl (boost::shared_ptr< ArcLengthCtx > &arc_ptr, const bool assemble=false)
 
 ~SimpleArcLengthControl ()
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode operator() ()
 
MoFEMErrorCode postProcess ()
 
double calculateLambdaInt ()
 Calculate internal lambda. More...
 
MoFEMErrorCode calculateDb ()
 Calculate db. More...
 
MoFEMErrorCode calculateDxAndDlambda (Vec x)
 

Public Attributes

boost::shared_ptr< ArcLengthCtxarcPtr
 
const bool aSsemble
 

Detailed Description

|brief Simple arc-length control of force

This is added for testing, it simply control force, i.e.

\[ \lambda = s \]

Constructor takes one argument,

Parameters
arc_ptrPointer to arc-length CTX.

Definition at line 334 of file ArcLengthTools.hpp.

Constructor & Destructor Documentation

◆ SimpleArcLengthControl()

SimpleArcLengthControl::SimpleArcLengthControl ( boost::shared_ptr< ArcLengthCtx > &  arc_ptr,
const bool  assemble = false 
)

Definition at line 521 of file ArcLengthTools.cpp.

523  : FEMethod(), arcPtr(arc_ptr), aSsemble(assemble) {}

◆ ~SimpleArcLengthControl()

SimpleArcLengthControl::~SimpleArcLengthControl ( )

Definition at line 525 of file ArcLengthTools.cpp.

525 {}

Member Function Documentation

◆ calculateDb()

MoFEMErrorCode SimpleArcLengthControl::calculateDb ( )

Calculate db.

Definition at line 597 of file ArcLengthTools.cpp.

597  {
599  CHKERR VecZeroEntries(arcPtr->db);
600  CHKERR VecGhostUpdateBegin(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
601  CHKERR VecGhostUpdateEnd(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
603 }

◆ calculateDxAndDlambda()

MoFEMErrorCode SimpleArcLengthControl::calculateDxAndDlambda ( Vec  x)

Definition at line 605 of file ArcLengthTools.cpp.

605  {
607  // Calculate dx
608  CHKERR VecGhostUpdateBegin(arcPtr->x0, INSERT_VALUES, SCATTER_FORWARD);
609  CHKERR VecGhostUpdateEnd(arcPtr->x0, INSERT_VALUES, SCATTER_FORWARD);
610 
611  Vec l_x, l_x0, l_dx;
612  CHKERR VecGhostGetLocalForm(x, &l_x);
613  CHKERR VecGhostGetLocalForm(arcPtr->x0, &l_x0);
614  CHKERR VecGhostGetLocalForm(arcPtr->dx, &l_dx);
615  {
616  double *x_array, *x0_array, *dx_array;
617  CHKERR VecGetArray(l_x, &x_array);
618  CHKERR VecGetArray(l_x0, &x0_array);
619  CHKERR VecGetArray(l_dx, &dx_array);
620  int size =
621  problemPtr->getNbLocalDofsRow() + problemPtr->getNbGhostDofsRow();
622  for (int i = 0; i != size; ++i) {
623  dx_array[i] = x_array[i] - x0_array[i];
624  }
625  CHKERR VecRestoreArray(l_x, &x_array);
626  CHKERR VecRestoreArray(l_x0, &x0_array);
627  CHKERR VecRestoreArray(l_dx, &dx_array);
628  }
629  CHKERR VecGhostRestoreLocalForm(x, &l_x);
630  CHKERR VecGhostRestoreLocalForm(arcPtr->x0, &l_x0);
631  CHKERR VecGhostRestoreLocalForm(arcPtr->dx, &l_dx);
632 
633  // Calculate dlambda
634  if (arcPtr->getPetscLocalDofIdx() != -1) {
635  double *array;
636  CHKERR VecGetArray(arcPtr->dx, &array);
637  arcPtr->dLambda = array[arcPtr->getPetscLocalDofIdx()];
638  array[arcPtr->getPetscLocalDofIdx()] = 0;
639  CHKERR VecRestoreArray(arcPtr->dx, &array);
640  }
641  CHKERR VecGhostUpdateBegin(arcPtr->ghosTdLambda, INSERT_VALUES,
642  SCATTER_FORWARD);
643  CHKERR VecGhostUpdateEnd(arcPtr->ghosTdLambda, INSERT_VALUES,
644  SCATTER_FORWARD);
645 
646  // Calculate dx2
647  double x_nrm, x0_nrm;
648  CHKERR VecNorm(x, NORM_2, &x_nrm);
649  CHKERR VecNorm(arcPtr->x0, NORM_2, &x0_nrm);
650  CHKERR VecDot(arcPtr->dx, arcPtr->dx, &arcPtr->dx2);
651 
652  MOFEM_LOG_C("WORLD", Sev::verbose,
653  "\tx norm = %6.4e x0 norm = %6.4e dx2 = %6.4e", x_nrm, x0_nrm,
654  arcPtr->dx2);
656 }

◆ calculateLambdaInt()

double SimpleArcLengthControl::calculateLambdaInt ( )

Calculate internal lambda.

Definition at line 593 of file ArcLengthTools.cpp.

593  {
594  return arcPtr->beta * arcPtr->dLambda;
595 }

◆ operator()()

MoFEMErrorCode SimpleArcLengthControl::operator() ( )

Definition at line 550 of file ArcLengthTools.cpp.

550  {
552  switch (snes_ctx) {
553  case CTX_SNESSETFUNCTION: {
554  arcPtr->res_lambda = calculateLambdaInt() - arcPtr->s;
555  CHKERR VecSetValue(snes_f, arcPtr->getPetscGlobalDofIdx(),
556  arcPtr->res_lambda, ADD_VALUES);
557  } break;
558  case CTX_SNESSETJACOBIAN: {
559  arcPtr->dIag = arcPtr->beta;
560  CHKERR MatSetValue(snes_B, arcPtr->getPetscGlobalDofIdx(),
561  arcPtr->getPetscGlobalDofIdx(), 1, ADD_VALUES);
562  } break;
563  default:
564  break;
565  }
567 }

◆ postProcess()

MoFEMErrorCode SimpleArcLengthControl::postProcess ( )

Definition at line 569 of file ArcLengthTools.cpp.

569  {
571  switch (snes_ctx) {
572  case CTX_SNESSETFUNCTION: {
573  if (aSsemble) {
574  CHKERR VecAssemblyBegin(snes_f);
575  CHKERR VecAssemblyEnd(snes_f);
576  }
577  } break;
578  case CTX_SNESSETJACOBIAN: {
579  if (aSsemble) {
580  CHKERR MatAssemblyBegin(snes_B, MAT_FLUSH_ASSEMBLY);
581  CHKERR MatAssemblyEnd(snes_B, MAT_FLUSH_ASSEMBLY);
582  }
583  CHKERR VecGhostUpdateBegin(arcPtr->ghostDiag, INSERT_VALUES,
584  SCATTER_FORWARD);
585  CHKERR VecGhostUpdateEnd(arcPtr->ghostDiag, INSERT_VALUES, SCATTER_FORWARD);
586  } break;
587  default:
588  break;
589  }
591 }

◆ preProcess()

MoFEMErrorCode SimpleArcLengthControl::preProcess ( )

Definition at line 527 of file ArcLengthTools.cpp.

527  {
529  switch (snes_ctx) {
530  case CTX_SNESSETFUNCTION: {
531  if (aSsemble) {
532  CHKERR VecAssemblyBegin(snes_f);
533  CHKERR VecAssemblyEnd(snes_f);
534  }
537  } break;
538  case CTX_SNESSETJACOBIAN: {
539  if (aSsemble) {
540  CHKERR MatAssemblyBegin(snes_B, MAT_FLUSH_ASSEMBLY);
541  CHKERR MatAssemblyEnd(snes_B, MAT_FLUSH_ASSEMBLY);
542  }
543  } break;
544  default:
545  break;
546  }
548 }

Member Data Documentation

◆ arcPtr

boost::shared_ptr<ArcLengthCtx> SimpleArcLengthControl::arcPtr

Definition at line 336 of file ArcLengthTools.hpp.

◆ aSsemble

const bool SimpleArcLengthControl::aSsemble

Definition at line 337 of file ArcLengthTools.hpp.


The documentation for this struct was generated from the following files:
ArcFunctionBegin
#define ArcFunctionBegin
Definition: ArcLengthTools.cpp:13
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
SimpleArcLengthControl::arcPtr
boost::shared_ptr< ArcLengthCtx > arcPtr
Definition: ArcLengthTools.hpp:336
SimpleArcLengthControl::aSsemble
const bool aSsemble
Definition: ArcLengthTools.hpp:337
SimpleArcLengthControl::calculateDxAndDlambda
MoFEMErrorCode calculateDxAndDlambda(Vec x)
Definition: ArcLengthTools.cpp:605
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
SimpleArcLengthControl::calculateDb
MoFEMErrorCode calculateDb()
Calculate db.
Definition: ArcLengthTools.cpp:597
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
SimpleArcLengthControl::calculateLambdaInt
double calculateLambdaInt()
Calculate internal lambda.
Definition: ArcLengthTools.cpp:593