Ice Sheet System Model  4.18
Code documentation
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Cfdragcoeffabsgrad.cpp
Go to the documentation of this file.
1 
5 /*Headers:*/
6 /*{{{*/
7 #ifdef HAVE_CONFIG_H
8  #include <config.h>
9 #else
10 #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
11 #endif
12 
13 #include "./classes.h"
16 #include "../datastructures/datastructures.h"
17 #include "./Elements/Element.h"
18 #include "./Elements/Elements.h"
19 #include "./FemModel.h"
20 #include "../modules/SurfaceAreax/SurfaceAreax.h"
21 #include "../classes/Params/Parameters.h"
22 #include "../classes/gauss/Gauss.h"
24 /*}}}*/
25 
26 /*Cfdragcoeffabsgrad constructors, destructors :*/
28 
29  this->definitionenum = -1;
30  this->name = NULL;
31  this->weights_enum = UNDEF;
32  this->misfit=0;
33  this->lock=0;
34  this->timepassedflag = false;
35  this->last_time = 0.;
36 
37 }
38 /*}}}*/
39 Cfdragcoeffabsgrad::Cfdragcoeffabsgrad(char* in_name, int in_definitionenum, int in_weights_enum, bool in_timepassedflag){/*{{{*/
40 
41  this->definitionenum=in_definitionenum;
42 
43  this->name = xNew<char>(strlen(in_name)+1);
44  xMemCpy<char>(this->name,in_name,strlen(in_name)+1);
45 
46  this->weights_enum=in_weights_enum;
47  this->timepassedflag=in_timepassedflag;
48 
49  this->misfit=0;
50  this->lock=0;
51 }
52 /*}}}*/
54  if(this->name)xDelete(this->name);
55  this->misfit=0;
56  this->lock=0;
57 }
58 /*}}}*/
59 /*Object virtual function resolutoin: */
62  mf->misfit=this->misfit;
63  mf->lock=this->lock;
64  return (Object*) mf;
65 }
66 /*}}}*/
67 void Cfdragcoeffabsgrad::DeepEcho(void){/*{{{*/
68  this->Echo();
69 }
70 /*}}}*/
71 void Cfdragcoeffabsgrad::Echo(void){/*{{{*/
72  _printf_(" Cfdragcoeffabsgrad: " << name << " " << this->definitionenum << "\n");
73  _printf_(" weights_enum: " << weights_enum << " " << EnumToStringx(weights_enum) << "\n");
74  _printf_(" timepassedflag: "<<timepassedflag<<"\n");
75 }
76 /*}}}*/
77 int Cfdragcoeffabsgrad::Id(void){/*{{{*/
78  return -1;
79 }
80 /*}}}*/
81 void Cfdragcoeffabsgrad::Marshall(char** pmarshalled_data,int* pmarshalled_data_size, int marshall_direction){/*{{{*/
82  _error_("not implemented yet!");
83 }
84 /*}}}*/
87 }
88 /*}}}*/
89 /*Definition virtual function resolutoin: */
91  return this->definitionenum;
92 }
93 /*}}}*/
94 char* Cfdragcoeffabsgrad::Name(){/*{{{*/
95  char* name2=xNew<char>(strlen(this->name)+1);
96  xMemCpy(name2,this->name,strlen(this->name)+1);
97 
98  return name2;
99 }
100 /*}}}*/
102  /*diverse: */
103  IssmDouble time;
104 
105  /*recover time parameters: */
106  int i;
107  IssmDouble J=0.;
108  IssmDouble J_sum=0.;
109 
110  for(i=0;i<femmodel->elements->Size();i++){
113  }
114 
115  ISSM_MPI_Allreduce ( (void*)&J,(void*)&J_sum,1,ISSM_MPI_DOUBLE,ISSM_MPI_SUM,IssmComm::GetComm());
117  J=J_sum;
118 
119  timepassedflag = true;
120  return J;
121  }
122  /*}}}*/
124 
125  int domaintype,numcomponents;
126  IssmDouble Jelem=0.;
127  IssmDouble misfit,Jdet;
128  IssmDouble dp[2],weight;
129  IssmDouble* xyz_list = NULL;
130 
131  /*Get basal element*/
132  if(!element->IsOnBase()) return 0.;
133 
134  /*If on water, return 0: */
135  if(!element->IsIceInElement()) return 0.;
136 
137  /*Get problem dimension*/
138  element->FindParam(&domaintype,DomainTypeEnum);
139  switch(domaintype){
140  case Domain2DverticalEnum: numcomponents = 1; break;
141  case Domain3DEnum: numcomponents = 2; break;
142  case Domain2DhorizontalEnum: numcomponents = 2; break;
143  default: _error_("not supported yet");
144  }
145 
146  /*Spawn surface element*/
147  Element* basalelement = element->SpawnBasalElement();
148 
149  /* Get node coordinates*/
150  basalelement->GetVerticesCoordinates(&xyz_list);
151 
152  /*Get input if it already exists*/
153  DatasetInput2 *datasetinput = basalelement->GetDatasetInput2(definitionenum); _assert_(datasetinput);
154  Input2 *drag_input = basalelement->GetInput2(FrictionCoefficientEnum); _assert_(drag_input);
155 
156  /* Start looping on the number of gaussian points: */
157  Gauss* gauss=basalelement->NewGauss(2);
158  for(int ig=gauss->begin();ig<gauss->end();ig++){
159 
160  gauss->GaussPoint(ig);
161 
162  /* Get Jacobian determinant: */
163  basalelement->JacobianDeterminant(&Jdet,xyz_list,gauss);
164 
165  /*Get all parameters at gaussian point*/
166  datasetinput->GetInputValue(&weight,gauss,WeightsSurfaceObservationEnum);
167  drag_input->GetInputDerivativeValue(&dp[0],xyz_list,gauss);
168 
169  /*Add to cost function*/
170  Jelem+=weight*.5*dp[0]*dp[0]*Jdet*gauss->weight;
171  if(numcomponents==2) Jelem+=weight*.5*dp[1]*dp[1]*Jdet*gauss->weight;
172  }
173 
174  /*clean up and Return: */
175  if(domaintype!=Domain2DhorizontalEnum){basalelement->DeleteMaterials(); delete basalelement;};
176  xDelete<IssmDouble>(xyz_list);
177  delete gauss;
178  return Jelem;
179 }/*}}}*/
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
Cfdragcoeffabsgrad::Cfdragcoeffabsgrad_Calculation
IssmDouble Cfdragcoeffabsgrad_Calculation(Element *element, int weights_enum)
Definition: Cfdragcoeffabsgrad.cpp:123
DatasetInput2
Definition: DatasetInput2.h:14
Cfdragcoeffabsgrad
Definition: Cfdragcoeffabsgrad.h:15
ISSM_MPI_Allreduce
int ISSM_MPI_Allreduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:116
Element::FindParam
void FindParam(bool *pvalue, int paramenum)
Definition: Element.cpp:933
_printf_
#define _printf_(StreamArgs)
Definition: Print.h:22
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
CfdragcoeffabsgradEnum
@ CfdragcoeffabsgradEnum
Definition: EnumDefinitions.h:1005
Element::SpawnBasalElement
virtual Element * SpawnBasalElement(void)=0
Cfdragcoeffabsgrad::DeepEcho
void DeepEcho(void)
Definition: Cfdragcoeffabsgrad.cpp:67
Cfdragcoeffabsgrad::weights_enum
int weights_enum
Definition: Cfdragcoeffabsgrad.h:21
Cfdragcoeffabsgrad::DefinitionEnum
int DefinitionEnum()
Definition: Cfdragcoeffabsgrad.cpp:90
Cfdragcoeffabsgrad::Echo
void Echo(void)
Definition: Cfdragcoeffabsgrad.cpp:71
DatasetInput2::GetInputValue
void GetInputValue(IssmDouble *pvalue, Gauss *gauss, int index)
Definition: DatasetInput2.cpp:199
Cfdragcoeffabsgrad::misfit
IssmDouble misfit
Definition: Cfdragcoeffabsgrad.h:26
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
Cfdragcoeffabsgrad::name
char * name
Definition: Cfdragcoeffabsgrad.h:20
Element
Definition: Element.h:41
Elements.h
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
Element.h
abstract class for Element object This class is a place holder for the Tria and the Penta elements....
Cfdragcoeffabsgrad::Cfdragcoeffabsgrad
Cfdragcoeffabsgrad()
Definition: Cfdragcoeffabsgrad.cpp:27
Object
Definition: Object.h:13
Cfdragcoeffabsgrad::timepassedflag
bool timepassedflag
Definition: Cfdragcoeffabsgrad.h:22
Cfdragcoeffabsgrad::copy
Object * copy()
Definition: Cfdragcoeffabsgrad.cpp:60
DomainTypeEnum
@ DomainTypeEnum
Definition: EnumDefinitions.h:124
xDelete
void xDelete(T **&aT_pp)
Definition: MemOps.h:97
Element::NewGauss
virtual Gauss * NewGauss(void)=0
Cfdragcoeffabsgrad::Response
IssmDouble Response(FemModel *femmodel)
Definition: Cfdragcoeffabsgrad.cpp:101
EnumToStringx
const char * EnumToStringx(int enum_in)
Definition: EnumToStringx.cpp:15
Cfdragcoeffabsgrad::Name
char * Name()
Definition: Cfdragcoeffabsgrad.cpp:94
ExternalResult.h
abstract class for ExternalResult object
UNDEF
#define UNDEF
Definition: constants.h:8
Element::GetVerticesCoordinates
void GetVerticesCoordinates(IssmDouble **xyz_list)
Definition: Element.cpp:1446
Cfdragcoeffabsgrad::last_time
IssmDouble last_time
Definition: Cfdragcoeffabsgrad.h:23
Domain3DEnum
@ Domain3DEnum
Definition: EnumDefinitions.h:536
FemModel::elements
Elements * elements
Definition: FemModel.h:44
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
FemModel
Definition: FemModel.h:31
Cfdragcoeffabsgrad::lock
int lock
Definition: Cfdragcoeffabsgrad.h:25
Cfdragcoeffabsgrad::definitionenum
int definitionenum
Definition: Cfdragcoeffabsgrad.h:19
Element::IsIceInElement
bool IsIceInElement()
Definition: Element.cpp:2021
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
Cfdragcoeffabsgrad::ObjectEnum
int ObjectEnum(void)
Definition: Cfdragcoeffabsgrad.cpp:85
Gauss::begin
virtual int begin(void)=0
DataSet::GetObjectByOffset
Object * GetObjectByOffset(int offset)
Definition: DataSet.cpp:334
FemModel.h
WeightsSurfaceObservationEnum
@ WeightsSurfaceObservationEnum
Definition: EnumDefinitions.h:864
Element::JacobianDeterminant
virtual void JacobianDeterminant(IssmDouble *Jdet, IssmDouble *xyz_list, Gauss *gauss)=0
Gauss::GaussPoint
virtual void GaussPoint(int ig)=0
Cfdragcoeffabsgrad::Marshall
void Marshall(char **pmarshalled_data, int *pmarshalled_data_size, int marshall_direction)
Definition: Cfdragcoeffabsgrad.cpp:81
Cfdragcoeffabsgrad::Id
int Id(void)
Definition: Cfdragcoeffabsgrad.cpp:77
Gauss::weight
IssmDouble weight
Definition: Gauss.h:11
xMemCpy
T * xMemCpy(T *dest, const T *src, unsigned int size)
Definition: MemOps.h:152
Domain2DverticalEnum
@ Domain2DverticalEnum
Definition: EnumDefinitions.h:535
Results.h
Gauss
Definition: Gauss.h:8
classes.h
DatasetInput2.h
: header file for datasetinput object
Cfdragcoeffabsgrad::~Cfdragcoeffabsgrad
~Cfdragcoeffabsgrad()
Definition: Cfdragcoeffabsgrad.cpp:53
femmodel
FemModel * femmodel
Definition: esmfbinders.cpp:16