Ice Sheet System Model  4.18
Code documentation
Functions
SurfaceAbsVelMisfitx.cpp File Reference
#include "./SurfaceAbsVelMisfitx.h"
#include "../../shared/shared.h"
#include "../../toolkits/toolkits.h"
#include "../../classes/Inputs2/DatasetInput2.h"

Go to the source code of this file.

Functions

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

Function Documentation

◆ SurfaceAbsVelMisfitx()

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

Definition at line 11 of file SurfaceAbsVelMisfitx.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+=SurfaceAbsVelMisfit(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 }

◆ SurfaceAbsVelMisfit()

IssmDouble SurfaceAbsVelMisfit ( Element element)

Definition at line 32 of file SurfaceAbsVelMisfitx.cpp.

32  {
33 
34  int domaintype,numcomponents;
35  IssmDouble Jelem=0.;
36  IssmDouble misfit,Jdet;
37  IssmDouble vx,vy,vxobs,vyobs,weight;
38  IssmDouble* xyz_list = NULL;
39 
40  /*Get basal element*/
41  if(!element->IsOnSurface()) 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 surface element*/
56  Element* topelement = element->SpawnTopElement();
57 
58  /* Get node coordinates*/
59  topelement->GetVerticesCoordinates(&xyz_list);
60 
61  /*Retrieve all inputs we will be needing: */
62  DatasetInput2* weights_input=topelement->GetDatasetInput2(InversionCostFunctionsCoefficientsEnum); _assert_(weights_input);
63  Input2* vx_input =topelement->GetInput2(VxEnum); _assert_(vx_input);
64  Input2* vxobs_input =topelement->GetInput2(InversionVxObsEnum); _assert_(vxobs_input);
65  Input2* vy_input = NULL;
66  Input2* vyobs_input = NULL;
67  if(numcomponents==2){
68  vy_input =topelement->GetInput2(VyEnum); _assert_(vy_input);
69  vyobs_input =topelement->GetInput2(InversionVyObsEnum); _assert_(vyobs_input);
70  }
71 
72  /* Start looping on the number of gaussian points: */
73  Gauss* gauss=topelement->NewGauss(2);
74  for(int ig=gauss->begin();ig<gauss->end();ig++){
75 
76  gauss->GaussPoint(ig);
77 
78  /* Get Jacobian determinant: */
79  topelement->JacobianDeterminant(&Jdet,xyz_list,gauss);
80 
81  /*Get all parameters at gaussian point*/
82  weights_input->GetInputValue(&weight,gauss,SurfaceAbsVelMisfitEnum);
83  vx_input->GetInputValue(&vx,gauss);
84  vxobs_input->GetInputValue(&vxobs,gauss);
85  if(numcomponents==2){
86  vy_input->GetInputValue(&vy,gauss);
87  vyobs_input->GetInputValue(&vyobs,gauss);
88  }
89 
90  /*Compute SurfaceAbsVelMisfitEnum:
91  *
92  * 1 [ 2 2 ]
93  * J = --- | (u - u ) + (v - v ) |
94  * 2 [ obs obs ]
95  *
96  */
97  misfit=0.5*(vx-vxobs)*(vx-vxobs);
98  if(numcomponents==2) misfit+=0.5*(vy-vyobs)*(vy-vyobs);
99 
100  /*Add to cost function*/
101  Jelem+=misfit*weight*Jdet*gauss->weight;
102  }
103 
104  /*clean up and Return: */
105  if(domaintype!=Domain2DhorizontalEnum){topelement->DeleteMaterials(); delete topelement;};
106  xDelete<IssmDouble>(xyz_list);
107  delete gauss;
108  return Jelem;
109 }
DataSet::Size
int Size()
Definition: DataSet.cpp:399
_assert_
#define _assert_(ignore)
Definition: exceptions.h:37
IssmDouble
double IssmDouble
Definition: types.h:37
DatasetInput2
Definition: DatasetInput2.h:14
Element::FindParam
void FindParam(bool *pvalue, int paramenum)
Definition: Element.cpp:933
InversionVyObsEnum
@ InversionVyObsEnum
Definition: EnumDefinitions.h:634
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
DatasetInput2::GetInputValue
void GetInputValue(IssmDouble *pvalue, Gauss *gauss, int index)
Definition: DatasetInput2.cpp:199
VyEnum
@ VyEnum
Definition: EnumDefinitions.h:850
SurfaceAbsVelMisfitEnum
@ SurfaceAbsVelMisfitEnum
Definition: EnumDefinitions.h:818
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
SurfaceAbsVelMisfit
IssmDouble SurfaceAbsVelMisfit(Element *element)
Definition: SurfaceAbsVelMisfitx.cpp:32
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
Element::SpawnTopElement
virtual Element * SpawnTopElement(void)=0
Gauss::begin
virtual int begin(void)=0
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
VxEnum
@ VxEnum
Definition: EnumDefinitions.h:846
Element::JacobianDeterminant
virtual void JacobianDeterminant(IssmDouble *Jdet, IssmDouble *xyz_list, Gauss *gauss)=0
Gauss::GaussPoint
virtual void GaussPoint(int ig)=0
Input2::GetInputValue
virtual void GetInputValue(IssmDouble *pvalue, Gauss *gauss)
Definition: Input2.h:38
InversionVxObsEnum
@ InversionVxObsEnum
Definition: EnumDefinitions.h:633
Gauss::weight
IssmDouble weight
Definition: Gauss.h:11
Element::IsOnSurface
bool IsOnSurface()
Definition: Element.cpp:1981
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