Changeset 16799
- Timestamp:
- 11/16/13 08:44:57 (11 years ago)
- Location:
- issm/trunk-jpl/src/c
- Files:
-
- 27 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/analyses/StressbalanceAnalysis.cpp
r16782 r16799 806 806 }/*}}}*/ 807 807 ElementVector* StressbalanceAnalysis::CreatePVector(Element* element){/*{{{*/ 808 _error_("not implemented yet"); 808 809 int approximation; 810 element->GetInputValue(&approximation,ApproximationEnum); 811 switch(approximation){ 812 case SSAApproximationEnum: 813 return CreatePVectorSSA(element); 814 default: 815 _error_("Approximation "<<EnumToStringx(approximation)<<" not supported"); 816 } 817 }/*}}}*/ 818 ElementVector* StressbalanceAnalysis::CreatePVectorSSA(Element* element){/*{{{*/ 819 820 /*compute all load vectors for this element*/ 821 ElementVector* pe1=CreatePVectorSSADrivingStress(element); 822 ElementVector* pe2=CreatePVectorSSAFront(element); 823 ElementVector* pe =new ElementVector(pe1,pe2); 824 825 /*clean-up and return*/ 826 delete pe1; 827 delete pe2; 828 return pe; 829 }/*}}}*/ 830 ElementVector* StressbalanceAnalysis::CreatePVectorSSADrivingStress(Element* element){/*{{{*/ 831 832 /*Intermediaries */ 833 int i; 834 IssmDouble driving_stress_baseline,thickness; 835 IssmDouble Jdet; 836 IssmDouble* xyz_list=NULL; 837 IssmDouble slope[2]; 838 839 /*Fetch number of nodes and dof for this finite element*/ 840 int numnodes = element->GetNumberOfNodes(); 841 IssmDouble* icefrontlevel= xNew<IssmDouble>(numnodes); 842 843 /*Initialize Element vector and vectors*/ 844 ElementVector* pe=element->NewElementVector(SSAApproximationEnum); 845 IssmDouble* basis = xNew<IssmDouble>(numnodes); 846 847 /*Retrieve all inputs and parameters*/ 848 element->GetVerticesCoordinates(&xyz_list); 849 Input* thickness_input=element->GetInput(ThicknessEnum); _assert_(thickness_input); 850 Input* surface_input =element->GetInput(SurfaceEnum); _assert_(surface_input); 851 Input* drag_input =element->GetInput(FrictionCoefficientEnum);_assert_(drag_input); 852 element->GetInputListOnVertices(icefrontlevel,BedEnum); 853 854 /* Start looping on the number of gaussian points: */ 855 Gauss* gauss=element->NewGauss(2); 856 for(int ig=gauss->begin();ig<gauss->end();ig++){ 857 858 gauss->GaussPoint(ig); 859 860 element->JacobianDeterminant(&Jdet,xyz_list,gauss); 861 element->NodalFunctions(basis, gauss); 862 863 thickness_input->GetInputValue(&thickness,gauss); 864 surface_input->GetInputDerivativeValue(&slope[0],xyz_list,gauss); 865 driving_stress_baseline=element->GetMaterialParameter(MaterialsRhoIceEnum)*element->GetMaterialParameter(ConstantsGEnum)*thickness; 866 867 /*Build load vector: */ 868 for(i=0;i<numnodes;i++){ 869 pe->values[i*NDOF2+0]+=-driving_stress_baseline*slope[0]*Jdet*gauss->weight*basis[i]; 870 pe->values[i*NDOF2+1]+=-driving_stress_baseline*slope[1]*Jdet*gauss->weight*basis[i]; 871 } 872 } 873 874 /*Transform coordinate system*/ 875 element->TransformLoadVectorCoord(pe,XYEnum); 876 877 /*Clean up and return*/ 878 xDelete<IssmDouble>(basis); 879 delete gauss; 880 return pe; 881 }/*}}}*/ 882 ElementVector* StressbalanceAnalysis::CreatePVectorSSAFront(Element* element){/*{{{*/ 883 809 884 }/*}}}*/ 810 885 void StressbalanceAnalysis::GetSolutionFromInputs(Vector<IssmDouble>* solution,Element* element){/*{{{*/ … … 952 1027 return; 953 1028 case L1L2ApproximationEnum: 954 InputUpdateFromSolution SSA(solution,element);1029 InputUpdateFromSolutionL1L2(solution,element); 955 1030 return; 956 1031 case SSAHOApproximationEnum: -
issm/trunk-jpl/src/c/analyses/StressbalanceAnalysis.h
r16782 r16799 23 23 ElementMatrix* CreateKMatrix(Element* element); 24 24 ElementVector* CreatePVector(Element* element); 25 ElementVector* CreatePVectorSSA(Element* element); 26 ElementVector* CreatePVectorSSADrivingStress(Element* element); 27 ElementVector* CreatePVectorSSAFront(Element* element); 25 28 void GetSolutionFromInputs(Vector<IssmDouble>* solution,Element* element); 26 29 void GetSolutionFromInputsFS(Vector<IssmDouble>* solution,Element* element); -
issm/trunk-jpl/src/c/classes/Elements/Element.h
r16789 r16799 24 24 class Input; 25 25 class Gauss; 26 class ElementVector; 26 27 template <class doublematrix> class Matrix; 27 28 template <class doubletype> class Vector; … … 54 55 virtual int FiniteElement(void)=0; 55 56 virtual IssmDouble GetMaterialParameter(int enum_in)=0; 57 virtual void NodalFunctions(IssmDouble* basis,Gauss* gauss)=0; 58 virtual void TransformLoadVectorCoord(ElementVector* pe,int transformenum)=0; 59 virtual void TransformLoadVectorCoord(ElementVector* pe,int* transformenum_list)=0; 60 virtual void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int transformenum)=0; 61 virtual void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int* transformenum_list)=0; 56 62 virtual void TransformSolutionCoord(IssmDouble* values,int transformenum)=0; 57 63 virtual void TransformSolutionCoord(IssmDouble* values,int* transformenum_list)=0; … … 62 68 virtual void GetDofListVelocity(int** pdoflist,int setenum)=0; 63 69 virtual void GetDofListPressure(int** pdoflist,int setenum)=0; 70 virtual void JacobianDeterminant(IssmDouble* Jdet, IssmDouble* xyz_list,Gauss* gauss)=0; 64 71 virtual void GetSolutionFromInputsOneDof(Vector<IssmDouble>* solution,int solutionenum)=0; 65 72 virtual int GetNodeIndex(Node* node)=0; … … 104 111 virtual int NumberofNodesPressure(void)=0; 105 112 virtual Gauss* NewGauss(void)=0; 113 virtual Gauss* NewGauss(int order)=0; 114 virtual ElementVector* NewElementVector(int approximation_enum)=0; 106 115 virtual void InputScale(int enum_type,IssmDouble scale_factor)=0; 107 116 virtual void GetVectorFromInputs(Vector<IssmDouble>* vector, int name_enum)=0; -
issm/trunk-jpl/src/c/classes/Elements/Penta.cpp
r16790 r16799 2519 2519 } 2520 2520 /*}}}*/ 2521 /*FUNCTION Penta::JacobianDeterminant{{{*/ 2522 void Penta::JacobianDeterminant(IssmDouble* pJdet,IssmDouble* xyz_list,Gauss* gauss){ 2523 2524 _assert_(gauss->Enum()==GaussPentaEnum); 2525 this->GetJacobianDeterminant(pJdet,xyz_list,(GaussPenta*)gauss); 2526 2527 } 2528 /*}}}*/ 2521 2529 /*FUNCTION Penta::NoIceInElement {{{*/ 2522 2530 bool Penta::NoIceInElement(){ … … 2587 2595 } 2588 2596 /*}}}*/ 2589 /*FUNCTION Penta::NewGauss {{{*/2597 /*FUNCTION Penta::NewGauss(){{{*/ 2590 2598 Gauss* Penta::NewGauss(void){ 2591 2599 return new GaussPenta(); 2600 } 2601 /*}}}*/ 2602 /*FUNCTION Penta::NewGauss(int order){{{*/ 2603 Gauss* Penta::NewGauss(int order){ 2604 return new GaussPenta(order,order); 2592 2605 } 2593 2606 /*}}}*/ … … 8021 8034 8022 8035 /*Transform coordinate system*/ 8023 TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list);8036 ::TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list); 8024 8037 8025 8038 /*Clean up and return*/ … … 8099 8112 8100 8113 /*Transform coordinate system*/ 8101 TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list);8114 ::TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list); 8102 8115 8103 8116 /*Clean up and return*/ … … 8183 8196 8184 8197 /*Transform coordinate system*/ 8185 TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list);8198 ::TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list); 8186 8199 8187 8200 /*Clean up and return*/ … … 8264 8277 8265 8278 /*Transform coordinate system*/ 8266 TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list);8279 ::TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list); 8267 8280 8268 8281 /*Clean up and return*/ … … 8542 8555 8543 8556 /*Transform coordinate system*/ 8544 TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum);8557 ::TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum); 8545 8558 8546 8559 /*Clean up and return*/ … … 8625 8638 8626 8639 /*Transform coordinate system*/ 8627 TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum);8640 ::TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum); 8628 8641 8629 8642 /*Clean up and return*/ … … 8733 8746 8734 8747 /*Transform coordinate system*/ 8735 TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list);8748 ::TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list); 8736 8749 8737 8750 /*Clean up and return*/ … … 8887 8900 8888 8901 /*Transform coordinate system*/ 8889 TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list);8902 ::TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list); 8890 8903 8891 8904 /*Clean up and return*/ … … 8972 8985 8973 8986 /*Transform coordinate system*/ 8974 TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list);8987 ::TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list); 8975 8988 8976 8989 /*Clean up and return*/ -
issm/trunk-jpl/src/c/classes/Elements/Penta.h
r16789 r16799 240 240 bool IsFloating(void); 241 241 bool IsNodeOnShelfFromFlags(IssmDouble* flags); 242 void JacobianDeterminant(IssmDouble* Jdet, IssmDouble* xyz_list,Gauss* gauss); 242 243 bool NoIceInElement(void); 243 244 Gauss* NewGauss(void); 245 Gauss* NewGauss(int order); 246 ElementVector* NewElementVector(int approximation_enum){_error_("not implemented yet");}; 247 void NodalFunctions(IssmDouble* basis,Gauss* gauss){_error_("not implemented yet");}; 244 248 IssmDouble MinEdgeLength(IssmDouble xyz_list[6][3]); 245 249 void SetClone(int* minranks); 246 250 Tria* SpawnTria(int location); 247 251 void SurfaceNormal(IssmDouble* surface_normal, IssmDouble xyz_list[3][3]); 252 void TransformLoadVectorCoord(ElementVector* pe,int transformenum){_error_("not implemented yet");}; 253 void TransformLoadVectorCoord(ElementVector* pe,int* transformenum_list){_error_("not implemented yet");}; 254 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int transformenum){_error_("not implemented yet");}; /*Tiling only*/ 255 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int* transformenum_list){_error_("not implemented yet");};/*Tiling only*/ 248 256 void TransformSolutionCoord(IssmDouble* values,int transformenum); 249 257 void TransformSolutionCoord(IssmDouble* values,int* transformenum_list); -
issm/trunk-jpl/src/c/classes/Elements/Seg.h
r16789 r16799 109 109 bool IsFloating(){_error_("not implemented yet");}; 110 110 bool IsNodeOnShelfFromFlags(IssmDouble* flags){_error_("not implemented yet");}; 111 void JacobianDeterminant(IssmDouble* Jdet, IssmDouble* xyz_list,Gauss* gauss){_error_("not implemented yet");}; 112 void NodalFunctions(IssmDouble* basis,Gauss* gauss){_error_("not implemented yet");}; 111 113 bool NoIceInElement(){_error_("not implemented yet");}; 112 114 int NumberofNodesVelocity(void){_error_("not implemented yet");}; … … 124 126 void GetMaterialInputValue(IssmDouble* pvalue,Node* node,int enumtype){_error_("not implemented yet");}; 125 127 Gauss* NewGauss(void){_error_("not implemented yet");}; 128 Gauss* NewGauss(int order){_error_("not implemented yet");}; 129 ElementVector* NewElementVector(int approximation_enum){_error_("not implemented yet");}; 126 130 #ifdef _HAVE_THERMAL_ 127 131 void UpdateBasalConstraintsEnthalpy(void){_error_("not implemented yet");}; … … 155 159 void Update(int index, IoModel* iomodel,int analysis_counter,int analysis_type,int finitelement){_error_("not implemented yet");}; 156 160 IssmDouble TimeAdapt(){_error_("not implemented yet");}; 161 void TransformLoadVectorCoord(ElementVector* pe,int transformenum){_error_("not implemented yet");}; 162 void TransformLoadVectorCoord(ElementVector* pe,int* transformenum_list){_error_("not implemented yet");}; 163 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int transformenum){_error_("not implemented yet");}; 164 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int* transformenum_list){_error_("not implemented yet");}; 157 165 void TransformSolutionCoord(IssmDouble* values,int transformenum){_error_("not implemented yet");}; 158 166 void TransformSolutionCoord(IssmDouble* values,int* transformenum_list){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Elements/Tria.cpp
r16789 r16799 2044 2044 } 2045 2045 /*}}}*/ 2046 /*FUNCTION Tria::JacobianDeterminant{{{*/ 2047 void Tria::JacobianDeterminant(IssmDouble* pJdet,IssmDouble* xyz_list,Gauss* gauss){ 2048 2049 _assert_(gauss->Enum()==GaussTriaEnum); 2050 this->JacobianDeterminant(pJdet,xyz_list,(GaussTria*)gauss); 2051 2052 } 2053 /*}}}*/ 2046 2054 /*FUNCTION Tria::HasEdgeOnBed {{{*/ 2047 2055 bool Tria::HasEdgeOnBed(){ … … 2229 2237 } 2230 2238 /*}}}*/ 2231 /*FUNCTION Tria::NewGauss {{{*/2239 /*FUNCTION Tria::NewGauss(){{{*/ 2232 2240 Gauss* Tria::NewGauss(void){ 2233 2241 return new GaussTria(); 2242 } 2243 /*}}}*/ 2244 /*FUNCTION Tria::NewGauss(int order){{{*/ 2245 Gauss* Tria::NewGauss(int order){ 2246 return new GaussTria(order); 2247 } 2248 /*}}}*/ 2249 /*FUNCTION Tria::NewElementVector{{{*/ 2250 ElementVector* Tria::NewElementVector(int approximation_enum){ 2251 return new ElementVector(nodes,this->NumberofNodes(),this->parameters,approximation_enum); 2234 2252 } 2235 2253 /*}}}*/ … … 2246 2264 /*If the level set is awlays <=0, there is no ice here*/ 2247 2265 return true; 2266 } 2267 /*}}}*/ 2268 /*FUNCTION Tria::NodalFunctions{{{*/ 2269 void Tria::NodalFunctions(IssmDouble* basis, Gauss* gauss){ 2270 2271 _assert_(gauss->Enum()==GaussTriaEnum); 2272 this->GetNodalFunctions(basis,(GaussTria*)gauss); 2273 2248 2274 } 2249 2275 /*}}}*/ … … 2731 2757 2732 2758 ::TransformSolutionCoord(values,this->nodes,this->NumberofNodes(),transformenum_list); 2759 2760 } 2761 /*}}}*/ 2762 /*FUNCTION Tria::TransformLoadVectorCoord{{{*/ 2763 void Tria::TransformLoadVectorCoord(ElementVector* pe,int transformenum){ 2764 2765 ::TransformLoadVectorCoord(pe,this->nodes,this->NumberofNodes(),transformenum); 2766 2767 } 2768 /*}}}*/ 2769 /*FUNCTION Tria::TransformLoadVectorCoord{{{*/ 2770 void Tria::TransformLoadVectorCoord(ElementVector* pe,int* transformenum_list){ 2771 2772 ::TransformLoadVectorCoord(pe,this->nodes,this->NumberofNodes(),transformenum_list); 2733 2773 2734 2774 } … … 3911 3951 3912 3952 /*Transform coordinate system*/ 3913 TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list);3953 ::TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list); 3914 3954 3915 3955 /*Clean up and return*/ … … 3972 4012 3973 4013 /*Transform coordinate system*/ 3974 TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list);4014 ::TransformLoadVectorCoord(pe,nodes,vnumnodes+pnumnodes,cs_list); 3975 4015 3976 4016 /*Clean up and return*/ … … 4042 4082 4043 4083 /*Transform coordinate system*/ 4044 TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list);4084 ::TransformLoadVectorCoord(pe,nodes,(vnumnodes+pnumnodes),cs_list); 4045 4085 4046 4086 /*Clean up and return*/ … … 4111 4151 4112 4152 /*Transform coordinate system*/ 4113 TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum);4153 ::TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum); 4114 4154 4115 4155 /*Clean up and return*/ … … 4191 4231 4192 4232 /*Transform coordinate system*/ 4193 TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum);4233 ::TransformLoadVectorCoord(pe,nodes,numnodes,XYEnum); 4194 4234 4195 4235 /*Clean up and return*/ -
issm/trunk-jpl/src/c/classes/Elements/Tria.h
r16789 r16799 268 268 void InputUpdateFromSolutionOneDofCollapsed(IssmDouble* solution,int enum_type){_error_("not implemented yet");}; 269 269 bool IsInput(int name); 270 void JacobianDeterminant(IssmDouble* pJdet, IssmDouble* xyz_list,Gauss* gauss); 270 271 Gauss* NewGauss(void); 272 Gauss* NewGauss(int order); 273 ElementVector* NewElementVector(int approximation_enum); 274 void NodalFunctions(IssmDouble* basis,Gauss* gauss); 271 275 void SetClone(int* minranks); 272 276 Seg* SpawnSeg(int index1,int index2); 273 277 void SurfaceNormal(IssmDouble* surface_normal, IssmDouble xyz_list[3][3]); 278 void TransformLoadVectorCoord(ElementVector* pe,int transformenum); 279 void TransformLoadVectorCoord(ElementVector* pe,int* transformenum_list); 280 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int transformenum){_error_("not implemented yet");}; 281 void TransformLoadVectorCoord(ElementVector* pe,int numnodes,int* transformenum_list){_error_("not implemented yet");}; 274 282 void TransformSolutionCoord(IssmDouble* values,int transformenum); 275 283 void TransformSolutionCoord(IssmDouble* values,int* transformenum_list); -
issm/trunk-jpl/src/c/classes/Inputs/BoolInput.cpp
r16675 r16799 115 115 void BoolInput::GetInputValue(IssmDouble* pvalue,Gauss* gauss){_error_("not supported yet!");} 116 116 /*}}}*/ 117 /*FUNCTION BoolInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){{{*/118 void BoolInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not supported yet!");}119 /*}}}*/120 /*FUNCTION BoolInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){{{*/121 void BoolInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not supported yet!");}122 /*}}}*/123 117 /*FUNCTION BoolInput::ChangeEnum{{{*/ 124 118 void BoolInput::ChangeEnum(int newenumtype){ -
issm/trunk-jpl/src/c/classes/Inputs/BoolInput.h
r16675 r16799 50 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 51 51 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 54 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss); 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, Gauss* gauss){_error_("not implemented yet");}; 55 53 void GetInputAverage(IssmDouble* pvalue){_error_("not implemented yet");}; 56 54 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/ControlInput.cpp
r16675 r16799 247 247 gradient->GetInputValue(pvalue,gauss); 248 248 }/*}}}*/ 249 /*FUNCTION ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){{{*/ 250 void ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){ 251 values->GetInputDerivativeValue(derivativevalues,xyz_list,gauss); 252 }/*}}}*/ 253 /*FUNCTION ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){{{*/ 254 void ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){ 249 /*FUNCTION ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, Gauss* gauss){{{*/ 250 void ControlInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, Gauss* gauss){ 255 251 values->GetInputDerivativeValue(derivativevalues,xyz_list,gauss); 256 252 }/*}}}*/ -
issm/trunk-jpl/src/c/classes/Inputs/ControlInput.h
r16675 r16799 56 56 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 57 57 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 58 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 59 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 60 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss); 58 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss); 61 59 void GetInputAverage(IssmDouble* pvalue); 62 60 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/DatasetInput.h
r16675 r16799 51 51 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 52 52 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index); 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 54 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not implemented yet");}; 55 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss){_error_("not implemented yet");}; 56 54 void GetInputAverage(IssmDouble* pvalue){_error_("not implemented yet");}; 57 55 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/DoubleInput.cpp
r16675 r16799 123 123 void DoubleInput::GetInputValue(IssmDouble* pvalue,Gauss* gauss){*pvalue=this->value;} 124 124 /*}}}*/ 125 /*FUNCTION DoubleInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){{{*/126 void DoubleInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not supported yet!");}127 /*}}}*/128 /*FUNCTION DoubleInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){{{*/129 void DoubleInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not supported yet!");}130 /*}}}*/131 125 /*FUNCTION DoubleInput::GetVxStrainRate2d(IssmDouble* epsilonvx,IssmDouble* xyz_list, GaussTria* gauss){{{*/ 132 126 void DoubleInput::GetVxStrainRate2d(IssmDouble* epsilonvx,IssmDouble* xyz_list, GaussTria* gauss){ -
issm/trunk-jpl/src/c/classes/Inputs/DoubleInput.h
r16675 r16799 49 49 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 50 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 51 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss); 51 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss){_error_("not implemented yet");}; 54 52 void GetInputAverage(IssmDouble* pvalue); 55 53 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/Input.h
r16675 r16799 32 32 virtual void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time)=0; 33 33 virtual void GetInputValue(IssmDouble* pvalue,Gauss* gauss,int index)=0; 34 virtual void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss)=0; 35 virtual void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss)=0; 36 virtual void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss)=0; 34 virtual void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, Gauss* gauss)=0; 37 35 virtual void GetInputAverage(IssmDouble* pvalue)=0; 38 36 virtual void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes)=0; -
issm/trunk-jpl/src/c/classes/Inputs/IntInput.cpp
r16675 r16799 115 115 void IntInput::GetInputValue(IssmDouble* pvalue,Gauss* gauss){_error_("not supported yet!");} 116 116 /*}}}*/ 117 /*FUNCTION IntInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){{{*/118 void IntInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not supported yet!");}119 /*}}}*/120 /*FUNCTION IntInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){{{*/121 void IntInput::GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not supported yet!");}122 /*}}}*/123 117 /*FUNCTION IntInput::ChangeEnum{{{*/ 124 118 void IntInput::ChangeEnum(int newenumtype){ -
issm/trunk-jpl/src/c/classes/Inputs/IntInput.h
r16675 r16799 50 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 51 51 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 54 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss); 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss){_error_("not implemented yet");}; 55 53 void GetInputAverage(IssmDouble* pvalue){_error_("not implemented yet");}; 56 54 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/PentaInput.cpp
r16675 r16799 146 146 } 147 147 /*}}}*/ 148 /*FUNCTION PentaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Penta* gauss){{{*/149 void PentaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, GaussPenta* gauss){148 /*FUNCTION PentaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){{{*/ 149 void PentaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list,Gauss* gauss){ 150 150 151 151 /*Call PentaRef function*/ 152 PentaRef::GetInputDerivativeValue(p,&values[0],xyz_list,gauss); 152 _assert_(gauss->Enum()==GaussPentaEnum); 153 PentaRef::GetInputDerivativeValue(p,&values[0],xyz_list,(GaussPenta*)gauss); 153 154 } 154 155 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Inputs/PentaInput.h
r16675 r16799 49 49 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 50 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 51 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss); 51 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss); 54 52 void GetInputAverage(IssmDouble* pvalue); 55 53 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/SegInput.cpp
r16675 r16799 105 105 } 106 106 /*}}}*/ 107 /*FUNCTION SegInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Seg* gauss){{{*/108 void SegInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, GaussSeg* gauss){107 /*FUNCTION SegInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){{{*/ 108 void SegInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list,Gauss* gauss){ 109 109 110 110 /*Call SegRef function*/ 111 SegRef::GetInputDerivativeValue(p,&values[0],xyz_list,gauss); 111 _assert_(gauss->Enum()==GaussSegEnum); 112 SegRef::GetInputDerivativeValue(p,&values[0],xyz_list,(GaussSeg*)gauss); 112 113 } 113 114 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Inputs/SegInput.h
r16675 r16799 48 48 void GetInputValue(IssmDouble* pvalue){_error_("not implemented yet");} 49 49 void GetInputValue(IssmDouble* pvalue,Gauss* gauss); 50 void GetInputValue(IssmDouble* pvalue,GaussTria* gauss){_error_("not implemented yet");};51 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 52 51 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss); 54 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss){_error_("not implemented yet");}; 55 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss); 56 53 void GetInputAverage(IssmDouble* pvalue); 57 54 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes){_error_("not implemented yet");}; -
issm/trunk-jpl/src/c/classes/Inputs/TransientInput.cpp
r16744 r16799 192 192 } 193 193 /*}}}*/ 194 /*FUNCTION TransientInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Tria* gauss){{{*/195 void TransientInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Tria* gauss){194 /*FUNCTION TransientInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){{{*/ 195 void TransientInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){ 196 196 197 197 IssmDouble time; … … 207 207 208 208 delete input; 209 210 209 } 211 210 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Inputs/TransientInput.h
r16740 r16799 55 55 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time); 56 56 void GetInputValue(IssmDouble* pvalue,Gauss* gauss ,int index){_error_("not implemented yet");}; 57 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 58 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 59 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not implemented yet");}; 57 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss); 60 58 void GetInputAverage(IssmDouble* pvalue); 61 59 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes); -
issm/trunk-jpl/src/c/classes/Inputs/TriaInput.cpp
r16675 r16799 133 133 } 134 134 /*}}}*/ 135 /*FUNCTION TriaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Tria* gauss){{{*/136 void TriaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss Tria* gauss){135 /*FUNCTION TriaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){{{*/ 136 void TriaInput::GetInputDerivativeValue(IssmDouble* p, IssmDouble* xyz_list, Gauss* gauss){ 137 137 138 138 /*Call TriaRef function*/ 139 TriaRef::GetInputDerivativeValue(p,&values[0],xyz_list,gauss); 139 _assert_(gauss->Enum()==GaussTriaEnum); 140 TriaRef::GetInputDerivativeValue(p,&values[0],xyz_list,(GaussTria*)gauss); 140 141 } 141 142 /*}}}*/ -
issm/trunk-jpl/src/c/classes/Inputs/TriaInput.h
r16675 r16799 50 50 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,IssmDouble time){_error_("not implemented yet");}; 51 51 void GetInputValue(IssmDouble* pvalue,Gauss* gauss,int index){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussSeg* gauss){_error_("not implemented yet");}; 53 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussTria* gauss); 54 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list, GaussPenta* gauss){_error_("not implemented yet");}; 52 void GetInputDerivativeValue(IssmDouble* derivativevalues, IssmDouble* xyz_list,Gauss* gauss); 55 53 void GetInputAverage(IssmDouble* pvalue); 56 54 void GetInputAllTimeAverages(IssmDouble** pvalues,IssmDouble** ptimes, int* pnumtimes); -
issm/trunk-jpl/src/c/classes/gauss/Gauss.h
r16675 r16799 9 9 10 10 public: 11 IssmDouble weight; 12 11 13 virtual ~Gauss(){}; 12 13 14 virtual int begin(void)=0; 14 15 virtual int end(void)=0;
Note:
See TracChangeset
for help on using the changeset viewer.