Ice Sheet System Model  4.18
Code documentation
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Functions
DragCoefficientAbsGradientx.h File Reference

header file for inverse methods misfit computation More...

#include "../../classes/classes.h"

Go to the source code of this file.

Functions

void DragCoefficientAbsGradientx (IssmDouble *pJ, Elements *elements, Nodes *nodes, Vertices *vertices, Loads *loads, Materials *materials, Parameters *parameters)
 
IssmDouble DragCoefficientAbsGradient (Element *element)
 

Detailed Description

header file for inverse methods misfit computation

Definition in file DragCoefficientAbsGradientx.h.

Function Documentation

◆ DragCoefficientAbsGradientx()

void DragCoefficientAbsGradientx ( IssmDouble pJ,
Elements elements,
Nodes nodes,
Vertices vertices,
Loads loads,
Materials materials,
Parameters parameters 
)

Definition at line 11 of file DragCoefficientAbsGradientx.cpp.

11  {
12 
13  /*output: */
14  IssmDouble J=0.;
15  IssmDouble J_sum;
16 
17  /*Compute Misfit: */
18  for(int i=0;i<elements->Size();i++){
19  Element* element=xDynamicCast<Element*>(elements->GetObjectByOffset(i));
20  J+=DragCoefficientAbsGradient(element);
21  }
22 
23  /*Sum all J from all cpus of the cluster:*/
26  J=J_sum;
27 
28  /*Assign output pointers: */
29  *pJ=J;
30 }

◆ DragCoefficientAbsGradient()

IssmDouble DragCoefficientAbsGradient ( Element element)

Definition at line 32 of file DragCoefficientAbsGradientx.cpp.

32  {
33 
34  int domaintype,numcomponents;
35  IssmDouble Jelem=0.;
36  IssmDouble misfit,Jdet;
37  IssmDouble dp[2],weight;
38  IssmDouble* xyz_list = NULL;
39 
40  /*Get basal element*/
41  if(!element->IsOnBase()) return 0.;
42 
43  /*If on water, return 0: */
44  if(!element->IsIceInElement()) return 0.;
45 
46  /*Get problem dimension*/
47  element->FindParam(&domaintype,DomainTypeEnum);
48  switch(domaintype){
49  case Domain2DverticalEnum: numcomponents = 1; break;
50  case Domain3DEnum: numcomponents = 2; break;
51  case Domain2DhorizontalEnum: numcomponents = 2; break;
52  default: _error_("not supported yet");
53  }
54 
55  /*Spawn basal element*/
56  Element* basalelement = element->SpawnBasalElement();
57 
58  /* Get node coordinates*/
59  basalelement->GetVerticesCoordinates(&xyz_list);
60 
61  /*Retrieve all inputs we will be needing: */
62  DatasetInput2* weights_input=basalelement->GetDatasetInput2(InversionCostFunctionsCoefficientsEnum); _assert_(weights_input);
63  Input2* drag_input =basalelement->GetInput2(FrictionCoefficientEnum); _assert_(drag_input);
64 
65  /* Start looping on the number of gaussian points: */
66  Gauss* gauss=basalelement->NewGauss(2);
67  for(int ig=gauss->begin();ig<gauss->end();ig++){
68 
69  gauss->GaussPoint(ig);
70 
71  /* Get Jacobian determinant: */
72  basalelement->JacobianDeterminant(&Jdet,xyz_list,gauss);
73 
74  /*Get all parameters at gaussian point*/
75  weights_input->GetInputValue(&weight,gauss,DragCoefficientAbsGradientEnum);
76  drag_input->GetInputDerivativeValue(&dp[0],xyz_list,gauss);
77 
78  /*Compute Tikhonov regularization J = 1/2 ((dp/dx)^2 + (dp/dy)^2) */
79  Jelem+=weight*.5*dp[0]*dp[0]*Jdet*gauss->weight;
80  if(numcomponents==2) Jelem+=weight*.5*dp[1]*dp[1]*Jdet*gauss->weight;
81 
82  }
83 
84  /*clean up and Return: */
85  if(domaintype!=Domain2DhorizontalEnum){basalelement->DeleteMaterials(); delete basalelement;};
86  xDelete<IssmDouble>(xyz_list);
87  delete gauss;
88  return Jelem;
89 }
DataSet::Size
int Size()
Definition: DataSet.cpp:399
FrictionCoefficientEnum
@ FrictionCoefficientEnum
Definition: EnumDefinitions.h:574
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IssmDouble
double IssmDouble
Definition: types.h:37
Element::IsOnBase
bool IsOnBase()
Definition: Element.cpp:1984
DatasetInput2
Definition: DatasetInput2.h:14
Element::FindParam
void FindParam(bool *pvalue, int paramenum)
Definition: Element.cpp:933
InversionCostFunctionsCoefficientsEnum
@ InversionCostFunctionsCoefficientsEnum
Definition: EnumDefinitions.h:629
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
Element::SpawnBasalElement
virtual Element * SpawnBasalElement(void)=0
DatasetInput2::GetInputValue
void GetInputValue(IssmDouble *pvalue, Gauss *gauss, int index)
Definition: DatasetInput2.cpp:199
DragCoefficientAbsGradientEnum
@ DragCoefficientAbsGradientEnum
Definition: EnumDefinitions.h:537
Input2::GetInputDerivativeValue
virtual void GetInputDerivativeValue(IssmDouble *derivativevalues, IssmDouble *xyz_list, Gauss *gauss)
Definition: Input2.h:37
Element::GetInput2
virtual Input2 * GetInput2(int inputenum)=0
Element::DeleteMaterials
void DeleteMaterials(void)
Definition: Element.cpp:429
Element
Definition: Element.h:41
Element::GetDatasetInput2
virtual DatasetInput2 * GetDatasetInput2(int inputenum)
Definition: Element.h:250
Domain2DhorizontalEnum
@ Domain2DhorizontalEnum
Definition: EnumDefinitions.h:534
ISSM_MPI_DOUBLE
#define ISSM_MPI_DOUBLE
Definition: issmmpi.h:125
DomainTypeEnum
@ DomainTypeEnum
Definition: EnumDefinitions.h:124
Element::NewGauss
virtual Gauss * NewGauss(void)=0
Element::GetVerticesCoordinates
void GetVerticesCoordinates(IssmDouble **xyz_list)
Definition: Element.cpp:1446
Domain3DEnum
@ Domain3DEnum
Definition: EnumDefinitions.h:536
ISSM_MPI_Bcast
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:162
Input2
Definition: Input2.h:18
Element::IsIceInElement
bool IsIceInElement()
Definition: Element.cpp:2021
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
Gauss::begin
virtual int begin(void)=0
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
Element::JacobianDeterminant
virtual void JacobianDeterminant(IssmDouble *Jdet, IssmDouble *xyz_list, Gauss *gauss)=0
Gauss::GaussPoint
virtual void GaussPoint(int ig)=0
DragCoefficientAbsGradient
IssmDouble DragCoefficientAbsGradient(Element *element)
Definition: DragCoefficientAbsGradientx.cpp:32
Gauss::weight
IssmDouble weight
Definition: Gauss.h:11
ISSM_MPI_Reduce
int ISSM_MPI_Reduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:373
Domain2DverticalEnum
@ Domain2DverticalEnum
Definition: EnumDefinitions.h:535
Gauss
Definition: Gauss.h:8