v0.14.0
Classes | Macros | Functions
TSEPAdapt.hpp File Reference

Step adaptation. More...

#include <petsc/private/tsimpl.h>

Go to the source code of this file.

Classes

struct  TSAdapt_EP
 

Macros

#define TSADAPTEP   "ep"
 

Functions

static PetscErrorCode TSAdaptChoose_EP (TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
 
static PetscErrorCode TSAdaptReset_EP (TSAdapt adapt)
 
static PetscErrorCode TSAdaptDestroy_EP (TSAdapt adapt)
 
PetscErrorCode TSAdaptCreate_EP (TSAdapt adapt)
 

Detailed Description

Step adaptation.

Problem implementation for mix element for large-strain elasticity

Todo:
Implementation of plasticity

Definition in file TSEPAdapt.hpp.

Macro Definition Documentation

◆ TSADAPTEP

#define TSADAPTEP   "ep"
Examples
ep.cpp.

Definition at line 15 of file TSEPAdapt.hpp.

Function Documentation

◆ TSAdaptChoose_EP()

static PetscErrorCode TSAdaptChoose_EP ( TSAdapt  adapt,
TS  ts,
PetscReal  h,
PetscInt *  next_sc,
PetscReal *  next_h,
PetscBool *  accept,
PetscReal *  wlte,
PetscReal *  wltea,
PetscReal *  wlter 
)
static

Definition at line 19 of file TSEPAdapt.hpp.

22  {
23  PetscFunctionBegin;
24 
25  TSAdapt_EP *basic = static_cast<TSAdapt_EP *>(adapt->data);
26 
27  *next_sc = 0; /* Reuse the same order scheme */
28  *wlte = -1; /* Weighted local truncation error was not evaluated */
29  *wltea = -1; /* Weighted absolute local truncation error is not used */
30  *wlter = -1; /* Weighted relative local truncation error is not used */
31 
32  *accept = PETSC_TRUE;
33  *next_h = h; /* Reuse the old step */
34 
35  SNES snes;
36  ierr = TSGetSNES(ts, &snes);
37  CHKERRG(ierr);
38 
39  SNESConvergedReason reason;
40  ierr = SNESGetConvergedReason(snes, &reason);
41  CHKERRG(ierr);
42 
43  int it;
44  ierr = SNESGetIterationNumber(snes, &it);
45  CHKERRG(ierr);
46 
47  if (reason < 0) {
48  h *= 0.75;
49  *next_h = h;
50  PetscPrintf(
51  PETSC_COMM_WORLD,
52  "\tDiverged set step length: it = %d, h = %3.4g set h = %3.4g \n", it,
53  h, *next_h);
54  } else if (reason > 0) {
55 
56  int its_d = 6;
57  ierr = PetscOptionsGetInt(PETSC_NULL, "", "-desired_nb_of_its", &its_d,
58  PETSC_NULL);
59  CHKERRG(ierr);
60 
61  h *= 1;//sqrt(static_cast<double>(its_d) / static_cast<double>(it + 1));
62  *next_h = PetscClipInterval(h, adapt->dt_min, adapt->dt_max);
63  PetscPrintf(
64  PETSC_COMM_WORLD,
65  "\tConverged set step length: it = %d, h = %3.4g set h = %3.4g \n", it,
66  h, *next_h);
67  }
68 
69  PetscFunctionReturn(0);
70 }

◆ TSAdaptCreate_EP()

PetscErrorCode TSAdaptCreate_EP ( TSAdapt  adapt)
Examples
ep.cpp.

Definition at line 87 of file TSEPAdapt.hpp.

87  {
88  TSAdapt_EP *ep;
89  PetscFunctionBegin;
90  ierr = PetscNewLog(adapt, &ep);
91  CHKERRG(ierr);
92  adapt->data = (void *)ep;
93  adapt->ops->choose = TSAdaptChoose_EP;
94  adapt->ops->reset = TSAdaptReset_EP;
95  adapt->ops->destroy = TSAdaptDestroy_EP;
96  PetscFunctionReturn(0);
97 }

◆ TSAdaptDestroy_EP()

static PetscErrorCode TSAdaptDestroy_EP ( TSAdapt  adapt)
static

Definition at line 78 of file TSEPAdapt.hpp.

78  {
79  PetscFunctionBegin;
80  ierr = TSAdaptReset_EP(adapt);
81  CHKERRG(ierr);
82  ierr = PetscFree(adapt->data);
83  CHKERRG(ierr);
84  PetscFunctionReturn(0);
85 }

◆ TSAdaptReset_EP()

static PetscErrorCode TSAdaptReset_EP ( TSAdapt  adapt)
static

Definition at line 72 of file TSEPAdapt.hpp.

72  {
73  TSAdapt_EP *basic = static_cast<TSAdapt_EP *>(adapt->data);
74  PetscFunctionBegin;
75  PetscFunctionReturn(0);
76 }
TSAdapt_EP
Definition: TSEPAdapt.hpp:17
TSAdaptReset_EP
static PetscErrorCode TSAdaptReset_EP(TSAdapt adapt)
Definition: TSEPAdapt.hpp:72
TSAdaptChoose_EP
static PetscErrorCode TSAdaptChoose_EP(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
Definition: TSEPAdapt.hpp:19
TSAdaptDestroy_EP
static PetscErrorCode TSAdaptDestroy_EP(TSAdapt adapt)
Definition: TSEPAdapt.hpp:78
h
double h
Definition: photon_diffusion.cpp:60
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496