v0.13.2
Loading...
Searching...
No Matches
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);
38
39 SNESConvergedReason reason;
40 ierr = SNESGetConvergedReason(snes,&reason);
42
43 int it;
44 ierr = SNESGetIterationNumber(snes, &it);
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 h *= sqrt(static_cast<double>(6) / static_cast<double>(it + 1));
56 *next_h = PetscClipInterval(h, adapt->dt_min, adapt->dt_max);
57 PetscPrintf(
58 PETSC_COMM_WORLD,
59 "\tConverged set step length: it = %d, h = %3.4g set h = %3.4g \n", it,
60 h, *next_h);
61 }
62
63 PetscFunctionReturn(0);
64}
static PetscErrorCode ierr
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
double h

◆ TSAdaptCreate_EP()

PetscErrorCode TSAdaptCreate_EP ( TSAdapt  adapt)
Examples
ep.cpp.

Definition at line 81 of file TSEPAdapt.hpp.

81 {
82 TSAdapt_EP *ep;
83 PetscFunctionBegin;
84 ierr = PetscNewLog(adapt, &ep);
86 adapt->data = (void *)ep;
87 adapt->ops->choose = TSAdaptChoose_EP;
88 adapt->ops->reset = TSAdaptReset_EP;
89 adapt->ops->destroy = TSAdaptDestroy_EP;
90 PetscFunctionReturn(0);
91}
static PetscErrorCode TSAdaptReset_EP(TSAdapt adapt)
Definition: TSEPAdapt.hpp:66
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
static PetscErrorCode TSAdaptDestroy_EP(TSAdapt adapt)
Definition: TSEPAdapt.hpp:72

◆ TSAdaptDestroy_EP()

static PetscErrorCode TSAdaptDestroy_EP ( TSAdapt  adapt)
static

Definition at line 72 of file TSEPAdapt.hpp.

72 {
73 PetscFunctionBegin;
74 ierr = TSAdaptReset_EP(adapt);
76 ierr = PetscFree(adapt->data);
78 PetscFunctionReturn(0);
79}

◆ TSAdaptReset_EP()

static PetscErrorCode TSAdaptReset_EP ( TSAdapt  adapt)
static

Definition at line 66 of file TSEPAdapt.hpp.

66 {
67 TSAdapt_EP *basic = static_cast<TSAdapt_EP *>(adapt->data);
68 PetscFunctionBegin;
69 PetscFunctionReturn(0);
70}